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