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.*
8 //
9 // -----------------------------
10 
12 
14 
16 
17 // -----------------------------
18 // constructors and destructor
19 // -----------------------------
20 
22  folderName_ ( iConfig.getParameter<std::string>("FolderName") )
23  , metToken_ ( consumes<reco::PFMETCollection> (iConfig.getParameter<edm::InputTag>("met") ) )
24  , jetToken_ ( mayConsume<reco::PFJetCollection> (iConfig.getParameter<edm::InputTag>("jets") ) )
25  , rsq_binning_ ( iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> > ("rsqBins") )
26  , mr_binning_ ( iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> > ("mrBins") )
27  , dphiR_binning_ ( iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> > ("dphiRBins") )
28  , num_genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"),consumesCollector(), *this))
29  , den_genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"),consumesCollector(), *this))
30  , metSelection_ ( iConfig.getParameter<std::string>("metSelection") )
31  , jetSelection_ ( iConfig.getParameter<std::string>("jetSelection") )
32  , njets_ ( iConfig.getParameter<unsigned int>("njets" ) )
33  , rsqCut_ ( iConfig.getParameter<double>("rsqCut" ) )
34  , mrCut_ ( iConfig.getParameter<double>("mrCut" ) )
35 {
36 
37  theHemispheres_ = consumes<std::vector<math::XYZTLorentzVector> >(iConfig.getParameter<edm::InputTag>("hemispheres"));
38 
39  MR_ME_.numerator = nullptr;
40  MR_ME_.denominator = nullptr;
41  Rsq_ME_.numerator = nullptr;
42  Rsq_ME_.denominator = nullptr;
43  dPhiR_ME_.numerator = nullptr;
44  dPhiR_ME_.denominator = nullptr;
45 
46  MRVsRsq_ME_.numerator = nullptr;
47  MRVsRsq_ME_.denominator = nullptr;
48 
49 }
50 
52 {
53 }
54 
56 {
57  me.numerator->setAxisTitle(titleX,1);
58  me.numerator->setAxisTitle(titleY,2);
59  me.denominator->setAxisTitle(titleX,1);
60  me.denominator->setAxisTitle(titleY,2);
61 
62 }
63 
64 void RazorMonitor::bookME(DQMStore::IBooker &ibooker, RazorME& me, const std::string& histname, const std::string& histtitle, int nbins, double min, double max)
65 {
66  me.numerator = ibooker.book1D(histname+"_numerator", histtitle+" (numerator)", nbins, min, max);
67  me.denominator = ibooker.book1D(histname+"_denominator", histtitle+" (denominator)", nbins, min, max);
68 }
69 void RazorMonitor::bookME(DQMStore::IBooker &ibooker, RazorME& me, const std::string& histname, const std::string& histtitle, const std::vector<double>& binning)
70 {
71  int nbins = binning.size()-1;
72  std::vector<float> fbinning(binning.begin(),binning.end());
73  float* arr = &fbinning[0];
74  me.numerator = ibooker.book1D(histname+"_numerator", histtitle+" (numerator)", nbins, arr);
75  me.denominator = ibooker.book1D(histname+"_denominator", histtitle+" (denominator)", nbins, arr);
76 }
77 void RazorMonitor::bookME(DQMStore::IBooker &ibooker, RazorME& me, const std::string& histname, const std::string& histtitle, int nbinsX, double xmin, double xmax, double ymin, double ymax)
78 {
79  me.numerator = ibooker.bookProfile(histname+"_numerator", histtitle+" (numerator)", nbinsX, xmin, xmax, ymin, ymax);
80  me.denominator = ibooker.bookProfile(histname+"_denominator", histtitle+" (denominator)", nbinsX, xmin, xmax, ymin, ymax);
81 }
82 void RazorMonitor::bookME(DQMStore::IBooker &ibooker, RazorME& me, const std::string& histname, const std::string& histtitle, int nbinsX, double xmin, double xmax, int nbinsY, double ymin, double ymax)
83 {
84  me.numerator = ibooker.book2D(histname+"_numerator", histtitle+" (numerator)", nbinsX, xmin, xmax, nbinsY, ymin, ymax);
85  me.denominator = ibooker.book2D(histname+"_denominator", histtitle+" (denominator)", nbinsX, xmin, xmax, nbinsY, ymin, ymax);
86 }
87 void RazorMonitor::bookME(DQMStore::IBooker &ibooker, RazorME& me, const std::string& histname, const std::string& histtitle, const std::vector<double>& binningX, const std::vector<double>& binningY)
88 {
89  int nbinsX = binningX.size()-1;
90  std::vector<float> fbinningX(binningX.begin(),binningX.end());
91  float* arrX = &fbinningX[0];
92  int nbinsY = binningY.size()-1;
93  std::vector<float> fbinningY(binningY.begin(),binningY.end());
94  float* arrY = &fbinningY[0];
95 
96  me.numerator = ibooker.book2D(histname+"_numerator", histtitle+" (numerator)", nbinsX, arrX, nbinsY, arrY);
97  me.denominator = ibooker.book2D(histname+"_denominator", histtitle+" (denominator)", nbinsX, arrX, nbinsY, arrY);
98 }
99 
101  edm::Run const & iRun,
102  edm::EventSetup const & iSetup)
103 {
104 
105  std::string histname, histtitle;
106 
107  std::string currentFolder = folderName_ ;
108  ibooker.setCurrentFolder(currentFolder.c_str());
109 
110  // 1D hist, MR
111  histname = "MR"; histtitle = "PF MR";
112  bookME(ibooker,MR_ME_,histname,histtitle,mr_binning_);
113  setMETitle(MR_ME_,"PF M_{R} [GeV]","events / [GeV]");
114 
115  // 1D hist, Rsq
116  histname = "Rsq"; histtitle = "PF Rsq";
117  bookME(ibooker,Rsq_ME_,histname,histtitle,rsq_binning_);
118  setMETitle(Rsq_ME_,"PF R^{2}","events");
119 
120  // 1D hist, dPhiR
121  histname = "dPhiR"; histtitle = "dPhiR";
122  bookME(ibooker,dPhiR_ME_,histname,histtitle,dphiR_binning_);
123  setMETitle(dPhiR_ME_,"dPhi_{R}","events");
124 
125  // 2D hist, MR & Rsq
126  histname = "MRVsRsq"; histtitle = "PF MR vs PF Rsq";
127  bookME(ibooker,MRVsRsq_ME_,histname,histtitle,mr_binning_, rsq_binning_);
128  setMETitle(MRVsRsq_ME_,"M_{R} [GeV]","R^{2}");
129 
130  // Initialize the GenericTriggerEventFlag
131  if ( num_genTriggerEventFlag_ && num_genTriggerEventFlag_->on() ) num_genTriggerEventFlag_->initRun( iRun, iSetup );
132  if ( den_genTriggerEventFlag_ && den_genTriggerEventFlag_->on() ) den_genTriggerEventFlag_->initRun( iRun, iSetup );
133 
134 }
135 
139 
140  // Filter out events if Trigger Filtering is requested
141  if (den_genTriggerEventFlag_->on() && ! den_genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
142 
143  //met collection
145  iEvent.getByToken( metToken_, metHandle );
146  reco::PFMET pfmet = metHandle->front();
147  if ( ! metSelection_( pfmet ) ) return;
148 
149  //jet collection, track # of jets for two working points
151  iEvent.getByToken( jetToken_, jetHandle );
152  std::vector<reco::PFJet> jets;
153  if ( jetHandle->size() < njets_ ) return;
154  for ( auto const & j : *jetHandle ) {
155  if ( jetSelection_( j ) ) jets.push_back(j);
156  }
157  if ( jets.size() < njets_ ) return;
158 
159  //razor hemisphere clustering from previous step
161  iEvent.getByToken (theHemispheres_,hemispheres);
162 
163  if (not hemispheres.isValid()){
164  return;
165  }
166 
167  if(hemispheres->size() ==0){ // the Hemisphere Maker will produce an empty collection of hemispheres if # of jets is too high
168  edm::LogError("DQM_HLT_Razor") << "Cannot calculate M_R and R^2 because there are too many jets! (trigger passed automatically without forming the hemispheres)" << endl;
169  return;
170  }
171 
172  // should always have 2 hemispheres -- no muons included (c. 2017), if not return invalid hemisphere collection
173  // retaining check for hemisphere size 5 or 10 which correspond to the one or two muon case for possible future use
174  if(hemispheres->size() != 0 && hemispheres->size() != 2 && hemispheres->size() != 5 && hemispheres->size() != 10){
175  edm::LogError("DQM_HLT_Razor") << "Invalid hemisphere collection! hemispheres->size() = " << hemispheres->size() << endl;
176  return;
177  }
178 
179  //calculate razor variables, with hemispheres pT-ordered
180  double MR = 0, R = 0;
181  if (hemispheres->at(1).Pt() > hemispheres->at(0).Pt()) {
182  MR = CalcMR(hemispheres->at(1),hemispheres->at(0));
183  R = CalcR(MR,hemispheres->at(1),hemispheres->at(0),metHandle);
184  }
185  else {
186  MR = CalcMR(hemispheres->at(0),hemispheres->at(1));
187  R = CalcR(MR,hemispheres->at(0),hemispheres->at(1),metHandle);
188  }
189 
190  double Rsq = R*R;
191  double dPhiR = abs(deltaPhi(hemispheres->at(0).Phi(), hemispheres->at(1).Phi()));
192 
193  //apply offline selection cuts
194  if (Rsq<rsqCut_ && MR<mrCut_) return;
195 
196  // applying selection for denominator
197  if (den_genTriggerEventFlag_->on() && ! den_genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
198 
199  // filling histograms (denominator)
200  if (Rsq>=rsqCut_) {
201  MR_ME_.denominator -> Fill(MR);
202  }
203 
204  if (MR>=mrCut_) {
205  Rsq_ME_.denominator -> Fill(Rsq);
206  }
207 
208  dPhiR_ME_.denominator -> Fill(dPhiR);
209 
210  MRVsRsq_ME_.denominator -> Fill(MR, Rsq);
211 
212  // applying selection for numerator
213  if (num_genTriggerEventFlag_->on() && ! num_genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
214 
215  // filling histograms (numerator)
216  if (Rsq>=rsqCut_) {
217  MR_ME_.numerator -> Fill(MR);
218  }
219 
220  if (MR>=mrCut_) {
221  Rsq_ME_.numerator -> Fill(Rsq);
222  }
223 
224  dPhiR_ME_.numerator -> Fill(dPhiR);
225 
226  MRVsRsq_ME_.numerator -> Fill(MR, Rsq);
227 
228 }
229 
231 {
233  desc.add<std::string> ( "FolderName", "HLT/SUSY/Razor" );
234 
235  desc.add<edm::InputTag>( "met", edm::InputTag("pfMet") );
236  desc.add<edm::InputTag>( "jets", edm::InputTag("ak4PFJetsCHS") );
237  desc.add<edm::InputTag>("hemispheres",edm::InputTag("hemispheresDQM"))->setComment("hemisphere jets used to compute razor variables");
238  desc.add<std::string>("metSelection", "pt > 0");
239 
240  // from 2016 offline selection
241  desc.add<std::string>("jetSelection", "pt > 80");
242  desc.add<unsigned int>("njets", 2);
243  desc.add<double>("mrCut", 300);
244  desc.add<double>("rsqCut", 0.15);
245 
246  edm::ParameterSetDescription genericTriggerEventPSet;
247  genericTriggerEventPSet.add<bool>("andOr");
248  genericTriggerEventPSet.add<edm::InputTag>("dcsInputTag", edm::InputTag("scalersRawToDigi") );
249  genericTriggerEventPSet.add<std::vector<int> >("dcsPartitions",{});
250  genericTriggerEventPSet.add<bool>("andOrDcs", false);
251  genericTriggerEventPSet.add<bool>("errorReplyDcs", true);
252  genericTriggerEventPSet.add<std::string>("dbLabel","");
253  genericTriggerEventPSet.add<bool>("andOrHlt", true);
254  genericTriggerEventPSet.add<edm::InputTag>("hltInputTag", edm::InputTag("TriggerResults::HLT") );
255  genericTriggerEventPSet.add<std::vector<std::string> >("hltPaths",{});
256  genericTriggerEventPSet.add<std::string>("hltDBKey","");
257  genericTriggerEventPSet.add<bool>("errorReplyHlt",false);
258  genericTriggerEventPSet.add<unsigned int>("verbosityLevel",1);
259 
260  desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
261  desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
262 
263  //binning from 2016 offline selection
265  std::vector<double> mrbins = {0., 100., 200., 300., 400., 500., 575., 650., 750., 900., 1200., 1600., 2500., 4000.};
266  histoPSet.add<std::vector<double> >("mrBins", mrbins);
267 
268  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};
269  histoPSet.add<std::vector<double> >("rsqBins", rsqbins);
270 
271  std::vector<double> dphirbins = {0., 0.5, 1.0, 1.5, 2.0, 2.5, 2.8, 3.0, 3.2};
272  histoPSet.add<std::vector<double> >("dphiRBins", dphirbins);
273 
274  desc.add<edm::ParameterSetDescription>("histoPSet",histoPSet);
275 
276  descriptions.add("razorMonitoring", desc);
277 }
278 
279 //CalcMR and CalcR borrowed from HLTRFilter.cc
281  if(ja.Pt()<=0.1) return -1;
282 
283  double A = ja.P();
284  double B = jb.P();
285  double az = ja.Pz();
286  double bz = jb.Pz();
287  TVector3 jaT, jbT;
288  jaT.SetXYZ(ja.Px(),ja.Py(),0.0);
289  jbT.SetXYZ(jb.Px(),jb.Py(),0.0);
290  double ATBT = (jaT+jbT).Mag2();
291 
292  double MR = sqrt((A+B)*(A+B)-(az+bz)*(az+bz)-
293  (jbT.Dot(jbT)-jaT.Dot(jaT))*(jbT.Dot(jbT)-jaT.Dot(jaT))/(jaT+jbT).Mag2());
294  double mybeta = (jbT.Dot(jbT)-jaT.Dot(jaT))/
295  sqrt(ATBT*((A+B)*(A+B)-(az+bz)*(az+bz)));
296 
297  double mygamma = 1./sqrt(1.-mybeta*mybeta);
298 
299  //use gamma times MRstar
300  return MR*mygamma;
301 }
302 
303 double RazorMonitor::CalcR(double MR, const math::XYZTLorentzVector& ja, const math::XYZTLorentzVector& jb, const edm::Handle<std::vector<reco::PFMET> >& inputMet){
304 
305  //now we can calculate MTR
306  const math::XYZVector met = (inputMet->front()).momentum();
307 
308  double MTR = sqrt(0.5*(met.R()*(ja.Pt()+jb.Pt()) - met.Dot(ja.Vect()+jb.Vect())));
309 
310  //filter events
311  return float(MTR)/float(MR); //R
312 }
313 
314 // Define this as a plug-in
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::PFJetCollection > jetToken_
Definition: RazorMonitor.h:85
edm::EDGetTokenT< reco::PFMETCollection > metToken_
Definition: RazorMonitor.h:84
void bookME(DQMStore::IBooker &, RazorME &me, const std::string &histname, const std::string &histtitle, int nbins, double xmin, double xmax)
Definition: RazorMonitor.cc:64
std::unique_ptr< GenericTriggerEventFlag > num_genTriggerEventFlag_
Definition: RazorMonitor.h:97
static double CalcMR(const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb)
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
MonitorElement * numerator
Definition: RazorMonitor.h:50
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
edm::EDGetTokenT< std::vector< math::XYZTLorentzVector > > theHemispheres_
Definition: RazorMonitor.h:86
unsigned int njets_
Definition: RazorMonitor.h:102
int iEvent
Definition: GenABIO.cc:230
std::vector< double > mr_binning_
Definition: RazorMonitor.h:89
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
static const std::string B
std::vector< double > dphiR_binning_
Definition: RazorMonitor.h:90
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
RazorME Rsq_ME_
Definition: RazorMonitor.h:93
void setMETitle(RazorME &me, std::string titleX, std::string titleY)
Definition: RazorMonitor.cc:55
RazorME dPhiR_ME_
Definition: RazorMonitor.h:94
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
std::vector< double > rsq_binning_
Definition: RazorMonitor.h:88
met
===> hadronic RAZOR
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
StringCutObjectSelector< reco::MET, true > metSelection_
Definition: RazorMonitor.h:100
RazorME MR_ME_
Definition: RazorMonitor.h:92
std::unique_ptr< GenericTriggerEventFlag > den_genTriggerEventFlag_
Definition: RazorMonitor.h:98
void add(std::string const &label, ParameterSetDescription const &psetDescription)
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
MonitorElement * denominator
Definition: RazorMonitor.h:51
RazorMonitor(const edm::ParameterSet &)
Definition: RazorMonitor.cc:21
static double CalcR(double MR, const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb, const edm::Handle< std::vector< reco::PFMET > > &met)
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
RazorME MRVsRsq_ME_
Definition: RazorMonitor.h:95
std::string folderName_
Definition: RazorMonitor.h:81
Definition: Run.h:42
StringCutObjectSelector< reco::PFJet, true > jetSelection_
Definition: RazorMonitor.h:101
Collection of PF MET.