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" ))
39 LogDebug(
"") <<
"Input/minJetPt/maxEta/maxNJ/acceptNJ : "
49 produces<std::vector<math::XYZTLorentzVector> >();
61 desc.
add<
bool>(
"doMuonCorrection",
false);
62 desc.
add<
double>(
"maxMuonEta",2.1);
63 desc.
add<
double>(
"minJetPt",30.0);
64 desc.
add<
double>(
"maxEta",3.0);
65 desc.
add<
int>(
"maxNJ",7);
66 desc.
add<
bool>(
"acceptNJ",
true);
67 descriptions.
add(
"hltRHemisphere",desc);
82 using namespace trigger;
96 std::auto_ptr<vector<math::XYZTLorentzVector> > Hemispheres(
new vector<math::XYZTLorentzVector> );
100 vector<math::XYZTLorentzVector> JETS;
101 for (
unsigned int i=0;
i<
jets->size();
i++) {
103 JETS.push_back(
jets->at(
i).p4());
109 iEvent.
put(Hemispheres);
115 int muonIndex[nMu] = { -1, -1 };
116 std::vector<reco::RecoChargedCandidate>::const_iterator muonIt;
119 for(muonIt =
muons->begin(); muonIt!=
muons->end(); muonIt++,index++){
122 iEvent.
put(Hemispheres);
125 muonIndex[nPassMu++] =
index;
131 std::vector<math::XYZTLorentzVector> muonJets;
133 muonJets.push_back(leadMu.
p4());
134 Hemispheres->push_back(leadMu.
p4());
139 muonJets.push_back(secondMu.
p4());
140 Hemispheres->push_back(secondMu.
p4());
142 muonJets.push_back(leadMu.
p4());
147 if(n<2)
return false;
154 iEvent.
put(Hemispheres);
160 std::vector<math::XYZTLorentzVector>* extraJets){
161 using namespace math;
162 using namespace reco;
165 int nJets = JETS.size();
166 if(extraJets) nJets+=extraJets->size();
169 hlist->push_back(j1R);
170 hlist->push_back(j2R);
173 unsigned int N_comb =
pow(2,nJets);
175 double M_minR = 9999999999.0;
176 unsigned int j_count;
177 for (
unsigned int i = 0;
i < N_comb;
i++) {
179 unsigned int itemp =
i;
181 unsigned int count = 0;
182 while (j_count > 0) {
183 if (itemp/j_count == 1){
184 if(count<JETS.size()) j_temp1 += JETS.at(count);
185 else j_temp1 +=extraJets->at(count-JETS.size());
187 if(count<JETS.size()) j_temp2 += JETS.at(count);
188 else j_temp2 +=extraJets->at(count-JETS.size());
190 itemp -= j_count * (itemp/j_count);
194 double M_temp = j_temp1.M2() + j_temp2.M2();
195 if (M_temp < M_minR) {
202 hlist->push_back(j1R);
203 hlist->push_back(j2R);
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
HLTRHemisphere(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< reco::RecoChargedCandidate > > m_theMuonToken
edm::EDGetTokenT< edm::View< reco::Jet > > m_theJetToken
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Abs< T >::type abs(const T &t)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void ComputeHemispheres(std::auto_ptr< std::vector< math::XYZTLorentzVector > > &hlist, const std::vector< math::XYZTLorentzVector > &JETS, std::vector< math::XYZTLorentzVector > *extraJets=0)
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Power< A, B >::type pow(const A &a, const B &b)
virtual bool filter(edm::Event &, const edm::EventSetup &)
math::PtEtaPhiELorentzVectorF LorentzVector