CMS 3D CMS Logo

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  produces<std::vector<double>>();
29 
30  // set Token(-s)
31  inputJetTagToken_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("inputJetTag"));
32 
33  LogDebug("") << "Inputs: " << inputJetTag_.encode() << " ";
34 }
35 
36 // ------------ method called to produce the data ------------
38  using namespace std;
39  using namespace edm;
40  using namespace reco;
41 
42  // get hold of collection of objects
44  iEvent.getByToken(inputJetTagToken_, calojet_handle);
45 
46  std::unique_ptr<std::vector<double>> result(new std::vector<double>);
47  // check the the input collections are available
48  if (calojet_handle.isValid()) {
49  std::vector<TLorentzVector> myJets;
50  reco::CaloJetCollection::const_iterator jetIt;
51  for (jetIt = calojet_handle->begin(); jetIt != calojet_handle->end(); jetIt++) {
52  TLorentzVector j;
53  j.SetPtEtaPhiE(jetIt->pt(), jetIt->eta(), jetIt->phi(), jetIt->energy());
54  myJets.push_back(j);
55  }
56 
57  double alphaT = CalcAlphaT(myJets);
58  double HT = CalcHT(myJets);
59 
60  result->push_back(alphaT);
61  result->push_back(HT);
62  }
63  iEvent.put(std::move(result));
64 }
65 
66 double AlphaTVarProducer::CalcAlphaT(const std::vector<TLorentzVector> &jets) const {
67  std::vector<double> ETs;
68  TVector3 MHT{CalcMHT(jets), 0.0, 0.0};
69  float HT = CalcHT(jets);
70  // float HT = 0;
71  for (unsigned int i = 0; i < jets.size(); i++) {
72  if (jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5)
73  ETs.push_back(jets[i].Et());
74  // HT += jets[i].Et();
75  }
76  if (ETs.size() < 2.)
77  return 0.0;
78  if (ETs.size() > 16.)
79  return 0.0;
80  float DHT = deltaHt(ETs);
81 
82  float AlphaT = alphaT(HT, DHT, MHT.Mag());
83 
84  return AlphaT;
85 }
86 
87 double AlphaTVarProducer::deltaHt(const std::vector<double> &ETs) {
88  if (ETs.size() > 16.)
89  return 9999999;
90  std::vector<double> diff(1 << (ETs.size() - 1), 0.);
91  for (unsigned i = 0; i < diff.size(); i++)
92  for (unsigned j = 0; j < ETs.size(); j++)
93  diff[i] += ETs[j] * (1 - 2 * (int(i >> j) & 1));
94  std::vector<double>::const_iterator it;
95  double min = 9999999;
96  for (it = diff.begin(); it != diff.end(); it++)
97  if (*it < min)
98  min = *it;
99  return min;
100 }
101 
102 double AlphaTVarProducer::alphaT(const double HT, const double DHT, const double MHT) {
103  return 0.5 * (HT - DHT) / sqrt(HT * HT - MHT * MHT);
104 }
105 
106 double AlphaTVarProducer::CalcHT(const std::vector<TLorentzVector> &jets) {
107  double HT = 0;
108  for (unsigned int i = 0; i < jets.size(); i++) {
109  if (jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5)
110  HT += jets[i].Et();
111  }
112 
113  return HT;
114 }
115 
116 double AlphaTVarProducer::CalcMHT(const std::vector<TLorentzVector> &jets) {
117  TVector3 MHT;
118  for (unsigned int i = 0; i < jets.size(); i++) {
119  if (jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5)
120  MHT -= jets[i].Vect();
121  }
122  MHT.SetZ(0.0);
123  return MHT.Mag();
124 }
125 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
static double CalcMHT(const std::vector< TLorentzVector > &)
std::string encode() const
Definition: InputTag.cc:159
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
static double alphaT(const double, const double, const double)
edm::InputTag inputJetTag_
AlphaTVarProducer(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::CaloJetCollection > inputJetTagToken_
static double deltaHt(const std::vector< double > &)
bool isValid() const
Definition: HandleBase.h:70
static double CalcHT(const std::vector< TLorentzVector > &)
fixed size matrix
HLT enums.
double CalcAlphaT(const std::vector< TLorentzVector > &) const
Definition: HT.h:21
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)