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.
3 
5 
7 
10 
13 
14 
17 
18 #include "TVector3.h"
19 #include "TLorentzVector.h"
20 
24 
25 #include<vector>
26 
27 //
28 // constructors and destructor
29 //
31  inputTag_ (iConfig.getParameter<edm::InputTag>("inputTag")),
32  min_Jet_Pt_ (iConfig.getParameter<double> ("minJetPt" )),
33  max_Eta_ (iConfig.getParameter<double> ("maxEta" )),
34  max_NJ_ (iConfig.getParameter<int> ("maxNJ" )),
35  accNJJets_ (iConfig.getParameter<bool> ("acceptNJ" ))
36 {
37  LogDebug("") << "Input/minJetPt/maxEta/maxNJ/acceptNJ : "
38  << inputTag_.encode() << " "
39  << min_Jet_Pt_ << "/"
40  << max_Eta_ << "/"
41  << max_NJ_ << "/"
42  << accNJJets_ << ".";
43 
44  //register your products
45  produces<std::vector<math::XYZTLorentzVector> >();
46 }
47 
49 {
50 }
51 
52 void
55  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltMCJetCorJetIcone5HF07"));
56  desc.add<double>("minJetPt",30.0);
57  desc.add<double>("maxEta",3.0);
58  desc.add<int>("maxNJ",7);
59  desc.add<bool>("acceptNJ",true);
60  descriptions.add("hltRHemisphere",desc);
61 }
62 
63 //
64 // member functions
65 //
66 
67 // ------------ method called to produce the data ------------
68 bool
70 {
71  using namespace std;
72  using namespace edm;
73  using namespace reco;
74  using namespace math;
75 
76  // get hold of collection of objects
78  iEvent.getByLabel (inputTag_,jets);
79 
80  // The output Collection
81  std::auto_ptr<vector<math::XYZTLorentzVector> > Hemispheres(new vector<math::XYZTLorentzVector> );
82 
83  // look at all objects, check cuts and add to filter object
84  int n(0);
86  CaloJetCollection::const_iterator i ( jets->begin() );
87  for (unsigned int i=0; i<jets->size(); i++) {
88  if(fabs(jets->at(i).eta()) < max_Eta_ && jets->at(i).pt() >= min_Jet_Pt_){
89  JETS.push_back(jets->at(i));
90  n++;
91  }
92  }
93 
94  if(n<2){
95  return false; //need at least 2 jets to build the hemispheres
96  }
97 
98  if(n>max_NJ_ && max_NJ_!=-1){
99  iEvent.put(Hemispheres);
100  return accNJJets_; //
101  }
102  int N_comb(1); // compute the number of combinations of jets possible
103  for(unsigned int i = 0; i < JETS.size(); i++){
104  N_comb *= 2;
105  }
106  //Make the hemispheres
107  XYZTLorentzVector j1,j2;
108  double M_min = 9999999999.0;
109  double dHT_min = 99999999.0;
110  int j_count;
111  for(int i=0;i<N_comb;i++){
112  XYZTLorentzVector j_temp1, j_temp2;
113  int itemp = i;
114  j_count = N_comb/2;
115  int count = 0;
116  while(j_count > 0){
117  if(itemp/j_count == 1){
118  j_temp1 += JETS.at(count).p4();
119  } else {
120  j_temp2 += JETS.at(count).p4();
121  }
122  itemp -= j_count*(itemp/j_count);
123  j_count /= 2;
124  count++;
125  }
126  double M_temp = j_temp1.M2()+j_temp2.M2();
127  if(M_temp < M_min){
128  M_min = M_temp;
129  j1= j_temp1;
130  j2= j_temp2;
131  }
132  double dHT_temp = fabs(j_temp1.E()-j_temp2.E());
133  if(dHT_temp < dHT_min){
134  dHT_min = dHT_temp;
135  }
136  }
137 
138  Hemispheres->push_back(j1);
139  Hemispheres->push_back(j2);
140 
141  iEvent.put(Hemispheres);
142  return true;
143 }
144 
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
HLTRHemisphere(const edm::ParameterSet &)
std::string encode() const
Definition: InputTag.cc:72
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:84
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::InputTag inputTag_
virtual bool filter(edm::Event &, const edm::EventSetup &)
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects