CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/JetMETCorrections/Modules/plugins/JetCorrectorDemo.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include "TH2F.h"
00003 #include "TGraph.h"
00004 #include "TGraphErrors.h"
00005 #include "TRandom.h"
00006 #include "TLorentzVector.h"
00007 // user include files
00008 #include "FWCore/Framework/interface/Frameworkfwd.h"
00009 #include "FWCore/Framework/interface/EDAnalyzer.h"
00010 
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/Framework/interface/MakerMacros.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "JetMETCorrections/Objects/interface/JetCorrector.h"
00017 #include "CondFormats/JetMETObjects/interface/JetCorrectionUncertainty.h"
00018 #include "JetMETCorrections/Objects/interface/JetCorrectionsRecord.h"
00019 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
00020 //TFile Service 
00021 #include "FWCore/ServiceRegistry/interface/Service.h"
00022 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00023 
00024 //
00025 // class declaration
00026 //
00027 
00028 class JetCorrectorDemo : public edm::EDAnalyzer {
00029 public:
00030   explicit JetCorrectorDemo(const edm::ParameterSet&);
00031   ~JetCorrectorDemo();
00032   typedef reco::Particle::LorentzVector LorentzVector;
00033   
00034 private:
00035   virtual void beginJob() ;
00036   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00037   virtual void endJob() ;
00038  
00039   std::string mJetCorService,mPayloadName,mUncertaintyTag,mUncertaintyFile;
00040   bool mDebug,mUseCondDB;
00041   int mNHistoPoints,mNGraphPoints;
00042   double mEtaMin,mEtaMax,mPtMin,mPtMax;
00043   std::vector<double> mVEta,mVPt;
00044   double vjec_eta[100][1000],vjec_pt[100][1000],vpt[100][1000],vptcor[100][1000],veta[100][1000];
00045   double vjecUnc_eta[100][1000],vUnc_eta[100][1000],vjecUnc_pt[100][1000],vUnc_pt[100][1000],vex_eta[100][1000],vex_pt[100][1000];
00046   edm::Service<TFileService> fs;
00047   TH2F *mJECvsEta, *mJECvsPt;
00048   TGraphErrors *mVGraphEta[100],*mVGraphPt[100],*mVGraphCorPt[100];
00049   TGraph *mUncEta[100], *mUncCorPt[100];
00050   TRandom *mRandom;
00051 };
00052 //
00053 //----------- Class Implementation ------------------------------------------
00054 //
00055 //---------------------------------------------------------------------------
00056 JetCorrectorDemo::JetCorrectorDemo(const edm::ParameterSet& iConfig)
00057 {
00058   mJetCorService     = iConfig.getParameter<std::string>          ("JetCorrectionService");
00059   mPayloadName       = iConfig.getParameter<std::string>          ("PayloadName");
00060   mUncertaintyTag    = iConfig.getParameter<std::string>          ("UncertaintyTag");
00061   mUncertaintyFile   = iConfig.getParameter<std::string>          ("UncertaintyFile");
00062   mNHistoPoints      = iConfig.getParameter<int>                  ("NHistoPoints");
00063   mNGraphPoints      = iConfig.getParameter<int>                  ("NGraphPoints");
00064   mEtaMin            = iConfig.getParameter<double>               ("EtaMin");
00065   mEtaMax            = iConfig.getParameter<double>               ("EtaMax");
00066   mPtMin             = iConfig.getParameter<double>               ("PtMin");
00067   mPtMax             = iConfig.getParameter<double>               ("PtMax");
00068   mVEta              = iConfig.getParameter<std::vector<double> > ("VEta");
00069   mVPt               = iConfig.getParameter<std::vector<double> > ("VPt");
00070   mDebug             = iConfig.getUntrackedParameter<bool>        ("Debug",false);
00071   mUseCondDB         = iConfig.getUntrackedParameter<bool>        ("UseCondDB",false);
00072 }
00073 //---------------------------------------------------------------------------
00074 JetCorrectorDemo::~JetCorrectorDemo()
00075 {
00076   
00077 }
00078 //---------------------------------------------------------------------------
00079 void JetCorrectorDemo::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00080 {
00081   const JetCorrector* corrector = JetCorrector::getJetCorrector(mJetCorService,iSetup);
00082   JetCorrectionUncertainty *jecUnc(0);
00083   if (mUncertaintyTag != "")
00084     {
00085       if (mUseCondDB)
00086         {
00087           edm::ESHandle<JetCorrectorParametersCollection> JetCorParColl;
00088           iSetup.get<JetCorrectionsRecord>().get(mPayloadName,JetCorParColl); 
00089           JetCorrectorParameters const & JetCorPar = (*JetCorParColl)[mUncertaintyTag];
00090           jecUnc = new JetCorrectionUncertainty(JetCorPar);
00091           std::cout<<"Configured Uncertainty from CondDB"<<std::endl;   
00092         }
00093       else
00094         {  
00095           edm::FileInPath fip("CondFormats/JetMETObjects/data/"+mUncertaintyFile+".txt");
00096           jecUnc = new JetCorrectionUncertainty(fip.fullPath()); 
00097         }
00098     }
00099   double jec,rawPt,corPt,eta,unc;
00100   TLorentzVector P4;
00101   double dEta = (mEtaMax-mEtaMin)/mNGraphPoints;
00102   for(int i=0;i<mNHistoPoints;i++)
00103     {
00104       rawPt  = mRandom->Uniform(mPtMin,mPtMax);
00105       eta = mRandom->Uniform(mEtaMin,mEtaMax);
00106       P4.SetPtEtaPhiE(rawPt,eta,0,0); 
00107       LorentzVector rawP4(P4.Px(),P4.Py(),P4.Pz(),P4.E());
00108       jec = corrector->correction(rawP4);
00109       mJECvsEta->Fill(eta,jec);
00110       mJECvsPt->Fill(rawPt,jec);
00111     }
00112   //--------- Pt Graphs ------------------
00113   for(unsigned ieta=0;ieta<mVEta.size();ieta++)
00114     {
00115       double rPt  = pow((3500./TMath::CosH(mVEta[ieta]))/mPtMin,1./mNGraphPoints);
00116       for(int i=0;i<mNGraphPoints;i++) 
00117         {
00118           rawPt  = mPtMin*pow(rPt,i);
00119           eta = mVEta[ieta];
00120           vpt[ieta][i] = rawPt;
00121           P4.SetPtEtaPhiE(rawPt,eta,0,0); 
00122           LorentzVector rawP4(P4.Px(),P4.Py(),P4.Pz(),P4.E());
00123           jec = corrector->correction(rawP4);// the jec is a function of the raw pt
00124           vjec_eta[ieta][i] = jec;
00125           vptcor[ieta][i] = rawPt*jec;
00126           unc = 0.0;
00127           if (mUncertaintyTag != "")
00128             {
00129               jecUnc->setJetEta(eta);
00130               jecUnc->setJetPt(rawPt*jec);// the uncertainty is a function of the corrected pt
00131               unc = jecUnc->getUncertainty(true);
00132             }
00133           vjecUnc_eta[ieta][i] = unc*jec;
00134           vUnc_eta[ieta][i] = unc;
00135           vex_eta[ieta][i] = 0.0;
00136           if (mDebug)
00137             std::cout<<rawPt<<" "<<eta<<" "<<jec<<" "<<rawPt*jec<<" "<<unc<<std::endl;
00138         }
00139     }
00140   //--------- Eta Graphs -------------------
00141   for(unsigned ipt=0;ipt<mVPt.size();ipt++)
00142     {
00143       for(int i=0;i<mNGraphPoints;i++)
00144         {
00145           eta = mEtaMin + i*dEta;
00146           corPt  = mVPt[ipt];
00147           veta[ipt][i] = eta;
00148           //---------- find the raw pt -----------
00149           double e = 1.0;
00150           int nLoop(0);
00151           rawPt = corPt; 
00152           while(e > 0.0001 && nLoop < 10) 
00153              {
00154                P4.SetPtEtaPhiE(rawPt,eta,0,0); 
00155                LorentzVector rawP4(P4.Px(),P4.Py(),P4.Pz(),P4.E()); 
00156                jec = corrector->correction(rawP4); 
00157                double tmp = rawPt * jec;
00158                e = fabs(tmp-corPt)/corPt;
00159                if (jec > 0)
00160                  rawPt = corPt/jec;
00161                nLoop++;
00162              } 
00163           //--------- calculate the jec for the rawPt --------
00164           P4.SetPtEtaPhiE(rawPt,eta,0,0); 
00165           LorentzVector rawP4(P4.Px(),P4.Py(),P4.Pz(),P4.E());
00166           jec = corrector->correction(rawP4);// the jec is a function of the raw pt
00167           vjec_pt[ipt][i] = jec;
00168           unc = 0.0;
00169           if (mUncertaintyTag != "")
00170             {
00171               jecUnc->setJetEta(eta);
00172               jecUnc->setJetPt(corPt);// the uncertainty is a function of the corrected pt
00173               unc = jecUnc->getUncertainty(true);
00174             }
00175           vjecUnc_pt[ipt][i] = unc*jec;
00176           vUnc_pt[ipt][i] = unc;
00177           vex_pt[ipt][i] = 0.0;
00178           if (mDebug)
00179             std::cout<<rawPt<<" "<<eta<<" "<<jec<<" "<<rawPt*jec<<" "<<unc<<std::endl;
00180         }
00181     }
00182 
00183 }
00184 //---------------------------------------------------------------------------
00185 void JetCorrectorDemo::beginJob()
00186 {
00187   if (mNGraphPoints > 1000)
00188     throw  cms::Exception("JetCorrectorDemo","Too many graph points !!! Maximum is 1000 !!!");
00189   if (mVEta.size() > 100)
00190     throw  cms::Exception("JetCorrectorDemo","Too many eta values !!! Maximum is 100 !!!");
00191   if (mVPt.size() > 100)
00192     throw  cms::Exception("JetCorrectorDemo","Too many pt values !!! Maximum is 100 !!!");
00193   mJECvsEta = fs->make<TH2F>("JECvsEta","JECvsEta",200,mEtaMin,mEtaMax,100,0,5);
00194   mJECvsPt  = fs->make<TH2F>("JECvsPt","JECvsPt",200,mPtMin,mPtMax,100,0,5);
00195   mRandom   = new TRandom();
00196   mRandom->SetSeed(0);
00197 }
00198 //---------------------------------------------------------------------------
00199 void JetCorrectorDemo::endJob() 
00200 {
00201   char name[1000];
00202   for(unsigned ipt=0;ipt<mVPt.size();ipt++)
00203     {
00204       mVGraphEta[ipt] = fs->make<TGraphErrors>(mNGraphPoints,veta[ipt],vjec_pt[ipt],vex_pt[ipt],vjecUnc_pt[ipt]);
00205       sprintf(name,"JEC_vs_Eta_CorPt%1.1f",mVPt[ipt]);
00206       mVGraphEta[ipt]->SetName(name);
00207       mUncEta[ipt] = fs->make<TGraph>(mNGraphPoints,veta[ipt],vUnc_pt[ipt]);
00208       sprintf(name,"UNC_vs_Eta_CorPt%1.1f",mVPt[ipt]);
00209       mUncEta[ipt]->SetName(name);
00210     }
00211   for(unsigned ieta=0;ieta<mVEta.size();ieta++)
00212     {
00213       mVGraphPt[ieta] = fs->make<TGraphErrors>(mNGraphPoints,vpt[ieta],vjec_eta[ieta],vex_eta[ieta],vjecUnc_eta[ieta]);
00214       sprintf(name,"JEC_vs_RawPt_eta%1.1f",mVEta[ieta]);
00215       mVGraphPt[ieta]->SetName(name);
00216       mVGraphCorPt[ieta] = fs->make<TGraphErrors>(mNGraphPoints,vptcor[ieta],vjec_eta[ieta],vex_eta[ieta],vjecUnc_eta[ieta]);
00217       sprintf(name,"JEC_vs_CorPt_eta%1.1f",mVEta[ieta]);
00218       mVGraphCorPt[ieta]->SetName(name);
00219       mUncCorPt[ieta] = fs->make<TGraph>(mNGraphPoints,vptcor[ieta],vUnc_eta[ieta]);
00220       sprintf(name,"UNC_vs_CorPt_eta%1.1f",mVEta[ieta]);
00221       mUncCorPt[ieta]->SetName(name);
00222     }
00223 }
00224 //---------------------------------------------------------------------------
00225 //define this as a plug-in
00226 DEFINE_FWK_MODULE(JetCorrectorDemo);
00227 
00228 
00229 
00230 
00231 
00232 
00233 
00234 
00235 
00236 
00237 
00238 
00239 
00240 
00241 
00242 
00243 
00244 
00245 
00246 
00247 
00248 
00249 
00250