41 static
double CalcR(
double MR,
91 requireValidHLTPaths_(
iConfig.getParameter<
bool>("requireValidHLTPaths")),
92 hltPathsAreValid_(
false),
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")),
108 rsqCut_(
iConfig.getParameter<
double>("rsqCut")),
109 mrCut_(
iConfig.getParameter<
double>("mrCut")) {}
151 histtitle =
"PF Rsq";
162 histname =
"MRVsRsq";
163 histtitle =
"PF MR vs PF Rsq";
189 std::vector<reco::PFJet>
jets;
190 if (jetHandle->size() <
njets_)
192 for (
auto const&
j : *jetHandle) {
203 if (not hemispheres.
isValid()) {
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)"
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()
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);
229 MR =
CalcMR(hemispheres->at(0), hemispheres->at(1));
230 R =
CalcR(MR, hemispheres->at(0), hemispheres->at(1), metHandle);
234 double dPhiR =
abs(
deltaPhi(hemispheres->at(0).Phi(), hemispheres->at(1).Phi()));
278 desc.
add<
bool>(
"requireValidHLTPaths",
true);
283 ->setComment(
"hemisphere jets used to compute razor variables");
288 desc.
add<
unsigned int>(
"njets", 2);
289 desc.
add<
double>(
"mrCut", 300);
290 desc.
add<
double>(
"rsqCut", 0.15);
293 genericTriggerEventPSet.
add<
bool>(
"andOr");
295 genericTriggerEventPSet.
add<std::vector<int> >(
"dcsPartitions", {});
296 genericTriggerEventPSet.
add<
bool>(
"andOrDcs",
false);
297 genericTriggerEventPSet.
add<
bool>(
"errorReplyDcs",
true);
299 genericTriggerEventPSet.
add<
bool>(
"andOrHlt",
true);
301 genericTriggerEventPSet.
add<std::vector<std::string> >(
"hltPaths", {});
303 genericTriggerEventPSet.
add<
bool>(
"errorReplyHlt",
false);
304 genericTriggerEventPSet.
add<
unsigned int>(
"verbosityLevel", 1);
311 std::vector<double> mrbins = {0., 100., 200., 300., 400., 500., 575., 650., 750., 900., 1200., 1600., 2500., 4000.};
312 histoPSet.
add<std::vector<double> >(
"mrBins", mrbins);
314 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};
315 histoPSet.add<std::vector<double> >(
"rsqBins", rsqbins);
317 std::vector<double> dphirbins = {0., 0.5, 1.0, 1.5, 2.0, 2.5, 2.8, 3.0, 3.2};
318 histoPSet.add<std::vector<double> >(
"dphiRBins", dphirbins);
322 descriptions.
add(
"razorMonitoring", desc);
335 jaT.SetXYZ(ja.Px(), ja.Py(), 0.0);
336 jbT.SetXYZ(jb.Px(), jb.Py(), 0.0);
337 double ATBT = (jaT + jbT).Mag2();
339 double MR =
sqrt((A + B) * (A + B) - (az + bz) * (az + bz) -
340 (jbT.Dot(jbT) - jaT.Dot(jaT)) * (jbT.Dot(jbT) - jaT.Dot(jaT)) / (jaT + jbT).Mag2());
341 double mybeta = (jbT.Dot(jbT) - jaT.Dot(jaT)) /
sqrt(ATBT * ((A + B) * (A + B) - (az + bz) * (az + bz)));
343 double mygamma = 1. /
sqrt(1. - mybeta * mybeta);
352 const edm::Handle<std::vector<reco::PFMET> >& inputMet) {
356 double MTR =
sqrt(0.5 * (met.R() * (ja.Pt() + jb.Pt()) - met.Dot(ja.Vect() + jb.Vect())));
359 return float(MTR) / float(MR);
dqm::reco::MonitorElement MonitorElement
edm::EDGetTokenT< reco::PFJetCollection > jetToken_
edm::EDGetTokenT< reco::PFMETCollection > metToken_
std::unique_ptr< GenericTriggerEventFlag > num_genTriggerEventFlag_
static double CalcMR(const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb)
virtual void setCurrentFolder(std::string const &fullpath)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
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)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
StringCutObjectSelector< reco::MET, true > metSelection_
edm::EDGetTokenT< std::vector< math::XYZTLorentzVector > > theHemispheres_
std::vector< double > mr_binning_
MonitorElement * denominator
MonitorElement * numerator
Abs< T >::type abs(const T &t)
static const std::string B
std::vector< double > dphiR_binning_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
StringCutObjectSelector< reco::PFJet, true > jetSelection_
std::vector< double > rsq_binning_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
std::unique_ptr< GenericTriggerEventFlag > den_genTriggerEventFlag_
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< PFJet > PFJetCollection
collection of PFJet objects
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
RazorMonitor(const edm::ParameterSet &)
static double CalcR(double MR, const math::XYZTLorentzVector &ja, const math::XYZTLorentzVector &jb, const edm::Handle< std::vector< reco::PFMET > > &met)
const std::string folderName_
const bool requireValidHLTPaths_
math::XYZTLorentzVector XYZTLorentzVector
Computes the MET from a collection of PFCandidates. HF missing!