CMS 3D CMS Logo

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