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  LogDebug("") << "Inputs: "
32  << inputJetTag_.encode() << " ";
33 }
34 
36 {
37 }
38 
39 // ------------ method called to produce the data ------------
40 void
42 {
43  using namespace std;
44  using namespace edm;
45  using namespace reco;
46 
47 
48  // get hold of collection of objects
50  iEvent.getByLabel(inputJetTag_,calojet_handle);
51 
52  std::auto_ptr<std::vector<double> > result(new std::vector<double>);
53  // check the the input collections are available
54  if (calojet_handle.isValid()){
55  std::vector<TLorentzVector> myJets;
56  reco::CaloJetCollection::const_iterator jetIt;
57  for(jetIt = calojet_handle->begin(); jetIt != calojet_handle->end(); jetIt++){
58  TLorentzVector j; j.SetPtEtaPhiE(jetIt->pt(),jetIt->eta(), jetIt->phi(), jetIt->energy());
59  myJets.push_back(j);
60  }
61 
62  double alphaT = CalcAlphaT(myJets);
63  double HT = CalcHT(myJets);
64 
65  result->push_back(alphaT);
66  result->push_back(HT);
67  }
68  iEvent.put(result);
69 }
70 
71 double
72 AlphaTVarProducer::CalcAlphaT(std::vector<TLorentzVector> jets){
73  std::vector<double> ETs;
74  TVector3 MHT = CalcMHT(jets);
75  float HT = CalcHT(jets);
76  //float HT = 0;
77  for(unsigned int i = 0; i < jets.size(); i++){
78  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) ETs.push_back(jets[i].Et());
79  //HT += jets[i].Et();
80  }
81  if(ETs.size() < 2.) return 0.0;
82  if(ETs.size() > 16.) return 0.0;
83  float DHT = deltaHt(ETs);
84 
85  float AlphaT = alphaT(HT,DHT,MHT.Mag());
86 
87  return AlphaT;
88 
89 }
90 
91 double AlphaTVarProducer::deltaHt(const std::vector<double>& ETs) {
92  if(ETs.size() > 16.) return 9999999;
93  std::vector<double> diff( 1<<(ETs.size()-1) , 0. );
94  for(unsigned i=0; i < diff.size(); i++)
95  for(unsigned j=0; j < ETs.size(); j++)
96  diff[i] += ETs[j] * ( 1 - 2 * (int(i>>j)&1) ) ;
97  std::vector<double>::const_iterator it;
98  double min=9999999;
99  for(it = diff.begin(); it !=diff.end(); it++) if(*it<min) min = *it;
100  return min;
101 }
102 
103 double AlphaTVarProducer::alphaT(const double HT, const double DHT, const double MHT) {
104  return 0.5 * ( HT - DHT ) / sqrt( HT*HT - MHT*MHT );
105 }
106 
107 double
108 AlphaTVarProducer::CalcHT(std::vector<TLorentzVector> jets){
109 
110  double HT = 0;
111  for(unsigned int i = 0; i < jets.size(); i++){
112  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) HT += jets[i].Et();
113  }
114 
115  return HT;
116 
117 }
118 
119 double AlphaTVarProducer::CalcMHT(std::vector<TLorentzVector> jets){
120  TVector3 MHT;
121  for(unsigned int i = 0; i < jets.size(); i++){
122  if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) MHT -= jets[i].Vect();
123  }
124  MHT.SetZ(0.0);
125  return MHT.Mag();
126 }
127 
#define LogDebug(id)
Definition: AlphaT.h:9
int i
Definition: DBlmapReader.cc:9
double CalcAlphaT(std::vector< TLorentzVector >)
static double CalcMHT(const std::vector< TLorentzVector >)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
static double alphaT(const double, const double, const double)
edm::InputTag inputJetTag_
#define min(a, b)
Definition: mlp_lapack.h:161
std::string encode() const
Definition: InputTag.cc:72
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:85
T sqrt(T t)
Definition: SSEVec.h:46
vector< PseudoJet > jets
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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