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