CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlphaTVarProducer.cc
Go to the documentation of this file.
2 
4 
7 
10 
13 
15 
17 
18 #include "TVector3.h"
19 
20 #include <memory>
21 #include <vector>
22 
23 //
24 // constructors and destructor
25 //
27  inputJetTag_ (iConfig.getParameter<edm::InputTag>("inputJetTag")){
28 
29  produces<std::vector<double> >();
30 
31  //set Token(-s)
32  inputJetTagToken_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("inputJetTag"));
33 
34  LogDebug("") << "Inputs: "
35  << inputJetTag_.encode() << " ";
36 }
37 
39 {
40 }
41 
42 // ------------ method called to produce the data ------------
43 void
45 {
46  using namespace std;
47  using namespace edm;
48  using namespace reco;
49 
50 
51  // get hold of collection of objects
53  iEvent.getByToken(inputJetTagToken_, calojet_handle);
54 
55  std::auto_ptr<std::vector<double> > result(new std::vector<double>);
56  // check the the input collections are available
57  if (calojet_handle.isValid()){
58  std::vector<TLorentzVector> myJets;
59  reco::CaloJetCollection::const_iterator jetIt;
60  for(jetIt = calojet_handle->begin(); jetIt != calojet_handle->end(); jetIt++){
61  TLorentzVector j; j.SetPtEtaPhiE(jetIt->pt(),jetIt->eta(), jetIt->phi(), jetIt->energy());
62  myJets.push_back(j);
63  }
64 
65  double alphaT = CalcAlphaT(myJets);
66  double HT = CalcHT(myJets);
67 
68  result->push_back(alphaT);
69  result->push_back(HT);
70  }
71  iEvent.put(result);
72 }
73 
74 double
75 AlphaTVarProducer::CalcAlphaT(const std::vector<TLorentzVector>& jets){
76  std::vector<double> ETs;
77  TVector3 MHT = CalcMHT(jets);
78  float HT = CalcHT(jets);
79  //float HT = 0;
80  for(unsigned int i = 0; i < jets.size(); i++){
81  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) ETs.push_back(jets[i].Et());
82  //HT += jets[i].Et();
83  }
84  if(ETs.size() < 2.) return 0.0;
85  if(ETs.size() > 16.) return 0.0;
86  float DHT = deltaHt(ETs);
87 
88  float AlphaT = alphaT(HT,DHT,MHT.Mag());
89 
90  return AlphaT;
91 
92 }
93 
94 double AlphaTVarProducer::deltaHt(const std::vector<double>& ETs) {
95  if(ETs.size() > 16.) return 9999999;
96  std::vector<double> diff( 1<<(ETs.size()-1) , 0. );
97  for(unsigned i=0; i < diff.size(); i++)
98  for(unsigned j=0; j < ETs.size(); j++)
99  diff[i] += ETs[j] * ( 1 - 2 * (int(i>>j)&1) ) ;
100  std::vector<double>::const_iterator it;
101  double min=9999999;
102  for(it = diff.begin(); it !=diff.end(); it++) if(*it<min) min = *it;
103  return min;
104 }
105 
106 double AlphaTVarProducer::alphaT(const double HT, const double DHT, const double MHT) {
107  return 0.5 * ( HT - DHT ) / sqrt( HT*HT - MHT*MHT );
108 }
109 
110 double
111 AlphaTVarProducer::CalcHT(const std::vector<TLorentzVector>& jets){
112 
113  double HT = 0;
114  for(unsigned int i = 0; i < jets.size(); i++){
115  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) HT += jets[i].Et();
116  }
117 
118  return HT;
119 
120 }
121 
122 double AlphaTVarProducer::CalcMHT(const std::vector<TLorentzVector>& jets){
123  TVector3 MHT;
124  for(unsigned int i = 0; i < jets.size(); i++){
125  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) MHT -= jets[i].Vect();
126  }
127  MHT.SetZ(0.0);
128  return MHT.Mag();
129 }
130 
#define LogDebug(id)
T getParameter(std::string const &) const
Definition: AlphaT.h:9
int i
Definition: DBlmapReader.cc:9
static double CalcMHT(const std::vector< TLorentzVector > &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static double alphaT(const double, const double, const double)
edm::InputTag inputJetTag_
double CalcAlphaT(const std::vector< TLorentzVector > &)
std::string encode() const
Definition: InputTag.cc:164
AlphaTVarProducer(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::CaloJetCollection > inputJetTagToken_
virtual void produce(edm::Event &, const edm::EventSetup &)
static double deltaHt(const std::vector< double > &)
static double CalcHT(const std::vector< TLorentzVector > &)
Definition: HT.h:20