CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTRHemisphere.cc
Go to the documentation of this file.
4 
7 
9 
12 
13 #include "TVector3.h"
14 #include "TLorentzVector.h"
15 
19 
21 
23 
24 #include<vector>
25 
26 //
27 // constructors and destructor
28 //
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" ))
38 {
39  LogDebug("") << "Input/minJetPt/maxEta/maxNJ/acceptNJ : "
40  << inputTag_.encode() << " "
41  << min_Jet_Pt_ << "/"
42  << max_Eta_ << "/"
43  << max_NJ_ << "/"
44  << accNJJets_ << ".";
45 
46  m_theJetToken = consumes<edm::View<reco::Jet>>(inputTag_);
47  m_theMuonToken = consumes<std::vector<reco::RecoChargedCandidate>>(muonTag_);
48  //register your products
49  produces<std::vector<math::XYZTLorentzVector> >();
50 }
51 
53 {
54 }
55 
56 void
59  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
60  desc.add<edm::InputTag>("muonTag",edm::InputTag(""));
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);
68 }
69 
70 //
71 // member functions
72 //
73 
74 // ------------ method called to produce the data ------------
75 bool
77 {
78  using namespace std;
79  using namespace edm;
80  using namespace reco;
81  using namespace math;
82  using namespace trigger;
83 
85 
86  // get hold of collection of objects
87  // Handle<CaloJetCollection> jets;
89  iEvent.getByToken (m_theJetToken,jets);
90 
91  // get hold of the muons, if necessary
94 
95  // The output Collection
96  std::auto_ptr<vector<math::XYZTLorentzVector> > Hemispheres(new vector<math::XYZTLorentzVector> );
97 
98  // look at all objects, check cuts and add to filter object
99  int n(0);
100  vector<math::XYZTLorentzVector> JETS;
101  for (unsigned int i=0; i<jets->size(); i++) {
102  if(std::abs(jets->at(i).eta()) < max_Eta_ && jets->at(i).pt() >= min_Jet_Pt_){
103  JETS.push_back(jets->at(i).p4());
104  n++;
105  }
106  }
107 
108  if(n>max_NJ_ && max_NJ_!=-1){
109  iEvent.put(Hemispheres);
110  return accNJJets_; // too many jets, accept for timing
111  }
112 
113  if(doMuonCorrection_){
114  const int nMu = 2;
115  int muonIndex[nMu] = { -1, -1 };
116  std::vector<reco::RecoChargedCandidate>::const_iterator muonIt;
117  int index = 0;
118  int nPassMu = 0;
119  for(muonIt = muons->begin(); muonIt!=muons->end(); muonIt++,index++){
120  if(std::abs(muonIt->eta()) > muonEta_ || muonIt->pt() < min_Jet_Pt_) continue; // skip muons out of eta range or too low pT
121  if(nPassMu >= 2){ // if we have already accepted two muons, accept the event
122  iEvent.put(Hemispheres); // too many muons, accept for timing
123  return true;
124  }
125  muonIndex[nPassMu++] = index;
126  }
127  //muons as MET
128  this->ComputeHemispheres(Hemispheres,JETS);
129  //lead muon as jet
130  if(nPassMu>0){
131  std::vector<math::XYZTLorentzVector> muonJets;
132  reco::RecoChargedCandidate leadMu = muons->at(muonIndex[0]);
133  muonJets.push_back(leadMu.p4());
134  Hemispheres->push_back(leadMu.p4());
135  this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // lead muon as jet
136  if(nPassMu>1){ // two passing muons
137  muonJets.pop_back();
138  reco::RecoChargedCandidate secondMu = muons->at(muonIndex[1]);
139  muonJets.push_back(secondMu.p4());
140  Hemispheres->push_back(secondMu.p4());
141  this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // lead muon as v, second muon as jet
142  muonJets.push_back(leadMu.p4());
143  this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // both muon as jets
144  }
145  }
146  }else{ // do MuonCorrection==false
147  if(n<2) return false; // not enough jets and not adding in muons
148  this->ComputeHemispheres(Hemispheres,JETS); // don't do the muon isolation, just run once and done
149  }
150  //Format:
151  // 0 muon: 2 hemispheres (2)
152  // 1 muon: 2 hemisheress + leadMuP4 + 2 hemispheres (5)
153  // 2 muon: 2 hemispheres + leadMuP4 + 2 hemispheres + 2ndMuP4 + 4 Hemispheres (10)
154  iEvent.put(Hemispheres);
155  return true;
156 }
157 
158 void
159 HLTRHemisphere::ComputeHemispheres(std::auto_ptr<std::vector<math::XYZTLorentzVector> >& hlist, const std::vector<math::XYZTLorentzVector>& JETS,
160  std::vector<math::XYZTLorentzVector>* extraJets){
161  using namespace math;
162  using namespace reco;
163  XYZTLorentzVector j1R(0.1, 0., 0., 0.1);
164  XYZTLorentzVector j2R(0.1, 0., 0., 0.1);
165  int nJets = JETS.size();
166  if(extraJets) nJets+=extraJets->size();
167 
168  if(nJets<2){ // put empty hemispheres if not enough jets
169  hlist->push_back(j1R);
170  hlist->push_back(j2R);
171  return;
172  }
173  unsigned int N_comb = pow(2,nJets); // compute the number of combinations of jets possible
174  //Make the hemispheres
175  double M_minR = 9999999999.0;
176  unsigned int j_count;
177  for (unsigned int i = 0; i < N_comb; i++) {
178  XYZTLorentzVector j_temp1, j_temp2;
179  unsigned int itemp = i;
180  j_count = N_comb/2;
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());
186  } else {
187  if(count<JETS.size()) j_temp2 += JETS.at(count);
188  else j_temp2 +=extraJets->at(count-JETS.size());
189  }
190  itemp -= j_count * (itemp/j_count);
191  j_count /= 2;
192  count++;
193  }
194  double M_temp = j_temp1.M2() + j_temp2.M2();
195  if (M_temp < M_minR) {
196  M_minR = M_temp;
197  j1R = j_temp1;
198  j2R = j_temp2;
199  }
200  }
201 
202  hlist->push_back(j1R);
203  hlist->push_back(j2R);
204  return;
205 }
206 
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
HLTRHemisphere(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< reco::RecoChargedCandidate > > m_theMuonToken
edm::EDGetTokenT< edm::View< reco::Jet > > m_theJetToken
std::string encode() const
Definition: InputTag.cc:164
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::InputTag muonTag_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
tuple muons
Definition: patZpeak.py:38
void ComputeHemispheres(std::auto_ptr< std::vector< math::XYZTLorentzVector > > &hlist, const std::vector< math::XYZTLorentzVector > &JETS, std::vector< math::XYZTLorentzVector > *extraJets=0)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
edm::InputTag inputTag_
virtual bool filter(edm::Event &, const edm::EventSetup &)
math::PtEtaPhiELorentzVectorF LorentzVector