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