14 #include "TLorentzVector.h" 30 : inputTag_(iConfig.getParameter<
edm::
InputTag>(
"inputTag")),
31 muonTag_(iConfig.getParameter<
edm::
InputTag>(
"muonTag")),
32 doMuonCorrection_(iConfig.getParameter<
bool>(
"doMuonCorrection")),
33 muonEta_(iConfig.getParameter<double>(
"maxMuonEta")),
34 min_Jet_Pt_(iConfig.getParameter<double>(
"minJetPt")),
35 max_Eta_(iConfig.getParameter<double>(
"maxEta")),
36 max_NJ_(iConfig.getParameter<
int>(
"maxNJ")),
37 accNJJets_(iConfig.getParameter<
bool>(
"acceptNJ")) {
44 produces<std::vector<math::XYZTLorentzVector>>();
53 desc.add<
bool>(
"doMuonCorrection",
false);
54 desc.add<
double>(
"maxMuonEta", 2.1);
55 desc.add<
double>(
"minJetPt", 30.0);
56 desc.add<
double>(
"maxEta", 3.0);
57 desc.add<
int>(
"maxNJ", 7);
58 desc.add<
bool>(
"acceptNJ",
true);
59 descriptions.
add(
"hltRHemisphere",
desc);
87 std::unique_ptr<vector<math::XYZTLorentzVector>> Hemispheres(
new vector<math::XYZTLorentzVector>);
91 vector<math::XYZTLorentzVector> JETS;
92 for (
auto const&
i : *
jets) {
94 JETS.push_back(
i.p4());
106 int muonIndex[nMu] = {-1, -1};
107 std::vector<reco::RecoChargedCandidate>::const_iterator muonIt;
110 for (muonIt =
muons->begin(); muonIt !=
muons->end(); muonIt++,
index++) {
117 muonIndex[nPassMu++] =
index;
123 std::vector<math::XYZTLorentzVector> muonJets;
125 muonJets.push_back(leadMu.
p4());
126 Hemispheres->push_back(leadMu.
p4());
131 muonJets.push_back(secondMu.
p4());
132 Hemispheres->push_back(secondMu.
p4());
134 muonJets.push_back(leadMu.
p4());
152 const std::vector<math::XYZTLorentzVector>& JETS,
153 std::vector<math::XYZTLorentzVector>* extraJets) {
154 using namespace math;
155 using namespace reco;
158 int nJets = JETS.size();
160 nJets += extraJets->size();
163 hlist->push_back(j1R);
164 hlist->push_back(j2R);
167 unsigned int N_comb =
pow(2, nJets);
169 double M_minR = 9999999999.0;
170 unsigned int j_count;
171 for (
unsigned int i = 0;
i < N_comb;
i++) {
173 unsigned int itemp =
i;
174 j_count = N_comb / 2;
175 unsigned int count = 0;
176 while (j_count > 0) {
177 if (itemp / j_count == 1) {
178 if (
count < JETS.size())
179 j_temp1 += JETS.at(
count);
181 j_temp1 += extraJets->at(
count - JETS.size());
183 if (
count < JETS.size())
184 j_temp2 += JETS.at(
count);
186 j_temp2 += extraJets->at(
count - JETS.size());
188 itemp -= j_count * (itemp / j_count);
192 double M_temp = j_temp1.M2() + j_temp2.M2();
193 if (M_temp < M_minR) {
200 hlist->push_back(j1R);
201 hlist->push_back(j2R);
#define DEFINE_FWK_MODULE(type)
HLTRHemisphere(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< reco::RecoChargedCandidate > > m_theMuonToken
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
edm::EDGetTokenT< edm::View< reco::Jet > > m_theJetToken
const LorentzVector & p4() const final
four-momentum Lorentz vector
void ComputeHemispheres(std::unique_ptr< std::vector< math::XYZTLorentzVector >> &hlist, const std::vector< math::XYZTLorentzVector > &JETS, std::vector< math::XYZTLorentzVector > *extraJets=nullptr)
Abs< T >::type abs(const T &t)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool filter(edm::Event &, const edm::EventSetup &) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Power< A, B >::type pow(const A &a, const B &b)
math::PtEtaPhiELorentzVectorF LorentzVector
math::XYZTLorentzVector XYZTLorentzVector
~HLTRHemisphere() override