CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/DQM/DataScouting/src/RazorVarProducer.cc

Go to the documentation of this file.
00001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00002 
00003 #include "DataFormats/Common/interface/Handle.h"
00004 
00005 #include "DataFormats/JetReco/interface/CaloJet.h"
00006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00007 
00008 #include "DataFormats/METReco/interface/CaloMET.h"
00009 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00010 
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 
00014 #include "FWCore/Utilities/interface/InputTag.h"
00015 
00016 #include "DQM/DataScouting/interface/RazorVarProducer.h"
00017 
00018 #include "TVector3.h"
00019 
00020 #include <memory>
00021 #include <vector>
00022 
00023 //
00024 // constructors and destructor
00025 //
00026 RazorVarProducer::RazorVarProducer(const edm::ParameterSet& iConfig) :
00027   inputTag_    (iConfig.getParameter<edm::InputTag>("inputTag")),
00028   inputMetTag_ (iConfig.getParameter<edm::InputTag>("inputMetTag")){
00029 
00030   produces<std::vector<double> >();
00031 
00032   LogDebug("") << "Inputs: "
00033                << inputTag_.encode() << " "
00034                << inputMetTag_.encode() << ".";
00035 }
00036 
00037 RazorVarProducer::~RazorVarProducer()
00038 {
00039 }
00040 
00041 // ------------ method called to produce the data  ------------ 
00042 void
00043 RazorVarProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00044 {
00045    using namespace std;
00046    using namespace edm;
00047    using namespace reco;
00048 
00049 
00050    // get hold of collection of objects
00051    Handle< vector<math::XYZTLorentzVector> > hemispheres;
00052    iEvent.getByLabel (inputTag_,hemispheres);
00053 
00054    // get hold of the MET Collection
00055    Handle<CaloMETCollection> inputMet;
00056    iEvent.getByLabel(inputMetTag_,inputMet);  
00057 
00058    std::auto_ptr<std::vector<double> > result(new std::vector<double>); 
00059    // check the the input collections are available
00060    if (hemispheres.isValid() && inputMet.isValid() && hemispheres->size() > 1){
00061 
00062      TLorentzVector ja(hemispheres->at(0).x(),hemispheres->at(0).y(),hemispheres->at(0).z(),hemispheres->at(0).t());
00063      TLorentzVector jb(hemispheres->at(1).x(),hemispheres->at(1).y(),hemispheres->at(1).z(),hemispheres->at(1).t());
00064 
00065      std::vector<math::XYZTLorentzVector> muonVec;
00066      const double MR = CalcMR(ja,jb);
00067      const double R  = CalcR(MR,ja,jb,inputMet,muonVec);
00068      result->push_back(MR);
00069      result->push_back(R);
00070      
00071    }
00072    iEvent.put(result);
00073 }
00074 
00075 double 
00076 RazorVarProducer::CalcMR(TLorentzVector ja, TLorentzVector jb){
00077   if(ja.Pt()<=0.1) return -1;
00078 
00079   ja.SetPtEtaPhiM(ja.Pt(),ja.Eta(),ja.Phi(),0.0);
00080   jb.SetPtEtaPhiM(jb.Pt(),jb.Eta(),jb.Phi(),0.0);
00081   
00082   if(ja.Pt() > jb.Pt()){
00083     TLorentzVector temp = ja;
00084     ja = jb;
00085     jb = temp;
00086   }
00087   
00088   double A = ja.P();
00089   double B = jb.P();
00090   double az = ja.Pz();
00091   double bz = jb.Pz();
00092   TVector3 jaT, jbT;
00093   jaT.SetXYZ(ja.Px(),ja.Py(),0.0);
00094   jbT.SetXYZ(jb.Px(),jb.Py(),0.0);
00095   double ATBT = (jaT+jbT).Mag2();
00096   
00097   double MR = sqrt((A+B)*(A+B)-(az+bz)*(az+bz)-
00098                    (jbT.Dot(jbT)-jaT.Dot(jaT))*(jbT.Dot(jbT)-jaT.Dot(jaT))/(jaT+jbT).Mag2());
00099   
00100   double mybeta = (jbT.Dot(jbT)-jaT.Dot(jaT))/
00101     sqrt(ATBT*((A+B)*(A+B)-(az+bz)*(az+bz)));
00102   
00103   double mygamma = 1./sqrt(1.-mybeta*mybeta);
00104   
00105   //use gamma times MRstar
00106   return MR*mygamma;  
00107 }
00108 
00109 double 
00110 RazorVarProducer::CalcR(double MR, TLorentzVector ja, TLorentzVector jb, edm::Handle<reco::CaloMETCollection> inputMet, std::vector<math::XYZTLorentzVector> muons){
00111   //now we can calculate MTR
00112   TVector3 met;
00113   met.SetPtEtaPhi((inputMet->front()).pt(),0.0,(inputMet->front()).phi());
00114   
00115   std::vector<math::XYZTLorentzVector>::const_iterator muonIt;
00116   for(muonIt = muons.begin(); muonIt!=muons.end(); muonIt++){
00117     TVector3 tmp;
00118     tmp.SetPtEtaPhi(muonIt->pt(),0,muonIt->phi());
00119     met-=tmp;
00120   }
00121 
00122   double MTR = sqrt(0.5*(met.Mag()*(ja.Pt()+jb.Pt()) - met.Dot(ja.Vect()+jb.Vect())));
00123   
00124   //filter events
00125   return float(MTR)/float(MR); //R
00126   
00127 }
00128 
00129 DEFINE_FWK_MODULE(RazorVarProducer);