CMS 3D CMS Logo

RazorVarProducer.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  : inputTag_(iConfig.getParameter<edm::InputTag>("inputTag")),
28  inputMetTag_(iConfig.getParameter<edm::InputTag>("inputMetTag")) {
29  // set Token(-s)
30  inputTagToken_ = consumes<std::vector<math::XYZTLorentzVector>>(iConfig.getParameter<edm::InputTag>("inputTag"));
31  inputMetTagToken_ = consumes<reco::CaloMETCollection>(iConfig.getParameter<edm::InputTag>("inputMetTag"));
32 
33  produces<std::vector<double>>();
34 
35  LogDebug("") << "Inputs: " << inputTag_.encode() << " " << inputMetTag_.encode() << ".";
36 }
37 
39 
40 // ------------ method called to produce the data ------------
42  using namespace std;
43  using namespace edm;
44  using namespace reco;
45 
46  // get hold of collection of objects
48  iEvent.getByToken(inputTagToken_, hemispheres);
49 
50  // get hold of the MET Collection
52  iEvent.getByToken(inputMetTagToken_, inputMet);
53 
54  std::unique_ptr<std::vector<double>> result(new std::vector<double>);
55  // check the the input collections are available
56  if (hemispheres.isValid() && inputMet.isValid() && hemispheres->size() > 1) {
57  TLorentzVector ja(hemispheres->at(0).x(), hemispheres->at(0).y(), hemispheres->at(0).z(), hemispheres->at(0).t());
58  TLorentzVector jb(hemispheres->at(1).x(), hemispheres->at(1).y(), hemispheres->at(1).z(), hemispheres->at(1).t());
59 
60  std::vector<math::XYZTLorentzVector> muonVec;
61  const double MR = CalcMR(ja, jb);
62  const double R = CalcR(MR, ja, jb, inputMet, muonVec);
63  result->push_back(MR);
64  result->push_back(R);
65  }
66  iEvent.put(std::move(result));
67 }
68 
69 double RazorVarProducer::CalcMR(TLorentzVector ja, TLorentzVector jb) {
70  if (ja.Pt() <= 0.1)
71  return -1;
72 
73  ja.SetPtEtaPhiM(ja.Pt(), ja.Eta(), ja.Phi(), 0.0);
74  jb.SetPtEtaPhiM(jb.Pt(), jb.Eta(), jb.Phi(), 0.0);
75 
76  if (ja.Pt() > jb.Pt()) {
77  TLorentzVector temp = ja;
78  ja = jb;
79  jb = temp;
80  }
81 
82  double A = ja.P();
83  double B = jb.P();
84  double az = ja.Pz();
85  double bz = jb.Pz();
86  TVector3 jaT, jbT;
87  jaT.SetXYZ(ja.Px(), ja.Py(), 0.0);
88  jbT.SetXYZ(jb.Px(), jb.Py(), 0.0);
89  double ATBT = (jaT + jbT).Mag2();
90 
91  double MR = sqrt((A + B) * (A + B) - (az + bz) * (az + bz) -
92  (jbT.Dot(jbT) - jaT.Dot(jaT)) * (jbT.Dot(jbT) - jaT.Dot(jaT)) / (jaT + jbT).Mag2());
93 
94  double mybeta = (jbT.Dot(jbT) - jaT.Dot(jaT)) / sqrt(ATBT * ((A + B) * (A + B) - (az + bz) * (az + bz)));
95 
96  double mygamma = 1. / sqrt(1. - mybeta * mybeta);
97 
98  // use gamma times MRstar
99  return MR * mygamma;
100 }
101 
102 double RazorVarProducer::CalcR(double MR,
103  const TLorentzVector &ja,
104  const TLorentzVector &jb,
106  const std::vector<math::XYZTLorentzVector> &muons) {
107  // now we can calculate MTR
108  TVector3 met;
109  met.SetPtEtaPhi((inputMet->front()).pt(), 0.0, (inputMet->front()).phi());
110 
111  std::vector<math::XYZTLorentzVector>::const_iterator muonIt;
112  for (muonIt = muons.begin(); muonIt != muons.end(); muonIt++) {
113  TVector3 tmp;
114  tmp.SetPtEtaPhi(muonIt->pt(), 0, muonIt->phi());
115  met -= tmp;
116  }
117 
118  double MTR = sqrt(0.5 * (met.Mag() * (ja.Pt() + jb.Pt()) - met.Dot(ja.Vect() + jb.Vect())));
119 
120  // filter events
121  return float(MTR) / float(MR); // R
122 }
123 
#define LogDebug(id)
T getParameter(std::string const &) const
RazorVarProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
double CalcMR(TLorentzVector ja, TLorentzVector jb)
std::string encode() const
Definition: InputTag.cc:159
edm::InputTag inputMetTag_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::CaloMETCollection > inputMetTagToken_
T sqrt(T t)
Definition: SSEVec.h:18
static const std::string B
bool isValid() const
Definition: HandleBase.h:74
met
===> hadronic RAZOR
double CalcR(double MR, const TLorentzVector &ja, const TLorentzVector &jb, edm::Handle< reco::CaloMETCollection > met, const std::vector< math::XYZTLorentzVector > &muons)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
fixed size matrix
HLT enums.
void produce(edm::Event &, const edm::EventSetup &) override
~RazorVarProducer() override
edm::InputTag inputTag_
edm::EDGetTokenT< std::vector< math::XYZTLorentzVector > > inputTagToken_
def move(src, dest)
Definition: eostools.py:511