CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Validation/RecoJets/plugins/CalibAnalyzer.h

Go to the documentation of this file.
00001 #ifndef CalibAnalyzer_h
00002 #define CalibAnalyzer_h
00003 
00004 #include <memory>
00005 #include <string>
00006 #include <vector>
00007 
00008 #include "TH1F.h"
00009 #include "TH2F.h"
00010 
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/Framework/interface/EDAnalyzer.h"
00013 #include "FWCore/Framework/interface/Frameworkfwd.h"
00014 #include "FWCore/Utilities/interface/InputTag.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 
00017 #include "FWCore/Utilities/interface/EDMException.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00020 
00021 #include "DataFormats/Math/interface/deltaR.h"
00022 #include "DataFormats/JetReco/interface/CaloJet.h"
00023 #include "Validation/RecoJets/interface/CompType.h"
00024 #include "Validation/RecoJets/interface/NameScheme.h"
00025 
00026 
00027 template <typename Ref, typename Rec, typename Alg>
00028 class CalibAnalyzer : public edm::EDAnalyzer {
00029 
00030  public:
00031 
00032   explicit CalibAnalyzer(const edm::ParameterSet&);
00033   ~CalibAnalyzer(){};
00034   
00035  private:
00036 
00037   virtual void beginJob();
00038   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00039   virtual void endJob(){ alg_.summarize(); };
00040 
00041   void fill(const double& var, const double& val, const std::vector<double>& bins, const std::vector<TH1F*>& hists);
00042 
00043  private:
00044 
00045   edm::InputTag recs_;
00046   edm::InputTag refs_;
00047   std::string hist_;
00048   int type_;                                      // comparison type (0:ratio/1:rel diff)
00049   int bins_;                                      // number of bins in fit histograms
00050   double min_, max_;                              // min/max of fit histograms
00051   std::vector<double> binsPt_, binsEta_;          // binning in pt/eta
00052 
00053  private:
00054   
00055   TH2F *recVsRef_;
00056   TH1F *calPt_, *resPt_;                          // cal/res vs pt
00057   std::vector<TH1F*> ktPt_;                       // calibration plots vs pt
00058   TH1F *calEta_, *resEta_;                        // cal/res vs eta
00059   std::vector<TH1F*> ktEta_;                      // calibration plots vs eta
00060   std::vector<TH1F*> calEtaPt_, resEtaPt_;        // cal/res vs eta and pt
00061   std::vector<std::vector<TH1F*> > ktEtaPt_;      // calibration plots vs eta and pt
00062 
00063   Alg alg_;                                       // matching/balancing algorithm
00064 };
00065 
00066 template <typename Ref, typename Rec, typename Alg>
00067 CalibAnalyzer<Ref, Rec, Alg>::CalibAnalyzer(const edm::ParameterSet& cfg):
00068   recs_( cfg.getParameter<edm::InputTag>( "recs" ) ),
00069   refs_( cfg.getParameter<edm::InputTag>( "refs" ) ),
00070   hist_( cfg.getParameter<std::string > ( "hist" ) ),
00071   type_( cfg.getParameter<int>( "type" ) ),
00072   bins_( cfg.getParameter<int>( "bins" ) ),
00073   min_ ( cfg.getParameter<double>( "min" ) ),
00074   max_ ( cfg.getParameter<double>( "max" ) ),
00075   binsPt_ ( cfg.getParameter<std::vector<double> >( "binsPt"  ) ),
00076   binsEta_( cfg.getParameter<std::vector<double> >( "binsEta" ) ),
00077   alg_( cfg )
00078 {
00079 }
00080 
00081 template <typename Ref, typename Rec, typename Alg>
00082 void CalibAnalyzer<Ref, Rec, Alg>::analyze(const edm::Event& evt, const edm::EventSetup& setup)
00083 {
00084   edm::Handle<Ref> refs;
00085   evt.getByLabel(refs_, refs);
00086 
00087   edm::Handle<Rec> recs;
00088   evt.getByLabel(recs_, recs);
00089   
00090   // do matching
00091   std::map<unsigned int, unsigned int> matches=alg_(*refs, *recs);
00092 
00093   if( !matches.size()>0 )
00094     edm::LogWarning ( "NoMatchOrBalance" ) 
00095       << "No single match/balance found to any Rec object in collection";
00096   
00097   // fill comparison plots for matched jets
00098   for(std::map<unsigned int, unsigned int>::const_iterator match=matches.begin(); 
00099       match!=matches.end(); ++match){
00100 
00101     CompType cmp(type_);
00102     double val=cmp((*refs)[match->first].pt(), (*recs)[match->second].pt());
00103 
00104     fill((*refs)[match->first].pt(),  val, binsPt_,  ktPt_ );// inclusive binned in pt 
00105     fill((*refs)[match->first].eta(), val, binsEta_, ktEta_);// inclusive binned in eta  
00106 
00107     // differential in eta binned in pt   
00108     for(int idx=0; idx<((int)binsEta_.size()-1); ++idx)
00109       if( (binsEta_[idx]<(*refs)[match->first].eta()) && ((*refs)[match->first].eta()<binsEta_[idx+1]) )
00110         fill((*refs)[match->first].pt(),  val, binsPt_,  ktEtaPt_[idx] );
00111     recVsRef_->Fill( TMath::Log10((*refs)[match->first].pt()), TMath::Log10((*recs)[match->second].pt()) );
00112   }
00113 }
00114 
00115 template <typename Ref, typename Rec, typename Alg>
00116 void CalibAnalyzer<Ref, Rec, Alg>::fill(const double& var, const double& val, const std::vector<double>& bins, const std::vector<TH1F*>& hists)
00117 {
00118   for(unsigned int idx=0; idx<(bins.size()-1); ++idx){
00119     if( (bins[idx]<var) && (var<bins[idx+1]) ){
00120       hists[idx]->Fill( val );
00121     }
00122   }
00123 }
00124 
00125 template <typename Ref, typename Rec, typename Alg>
00126 void CalibAnalyzer<Ref, Rec, Alg>::beginJob()
00127 {
00128   if( hist_.empty() )
00129     return;
00130 
00131   edm::Service<TFileService> fs;
00132   if( !fs )
00133     throw edm::Exception( edm::errors::Configuration, "TFile Service is not registered in cfg file" );
00134 
00135   ofstream hist(hist_.c_str(), std::ios::out);
00136   NameScheme val("val"), fit("fit"), cal("cal"), res("res");
00137 
00138   // book additional control histograms
00139   recVsRef_= fs->make<TH2F>( val.name("recVsRef"),val.name("recVsRef"), 20, 1., 3., 20, 1., 3.);
00140 
00141   // book kt histograms differential in pt
00142   for(int idx=0; idx<((int)binsPt_.size()-1); ++idx)
00143     ktPt_.push_back( fs->make<TH1F>(fit.name(hist, "ktPt",idx), fit.name("kt",idx), bins_, min_, max_) );
00144   calPt_= fs->make<TH1F>(cal.name(hist, "ktPt"), cal.name("calPt"), ((int)binsPt_.size()-1), &binsPt_[0]);
00145   resPt_= fs->make<TH1F>(res.name(hist, "ktPt"), res.name("resPt"), ((int)binsPt_.size()-1), &binsPt_[0]);
00146   
00147   // book kt histograms differential in eta  
00148   for(int jdx=0; jdx<((int)binsEta_.size()-1); ++jdx)
00149     ktEta_.push_back(fs->make<TH1F>(fit.name(hist, "ktEta",jdx),fit.name("kt",jdx), bins_, min_, max_) );
00150   calEta_= fs->make<TH1F>(cal.name(hist, "ktEta"), cal.name("calEta"), ((int)binsEta_.size()-1), &binsEta_[0]);
00151   resEta_= fs->make<TH1F>(res.name(hist, "ktEta"), res.name("resEta"), ((int)binsEta_.size()-1), &binsEta_[0]);
00152   
00153   // book kt histograms differential in eta and pt
00154   for(int jdx=0; jdx<((int)binsEta_.size()-1); ++jdx){
00155     std::vector<TH1F*> buffer;
00156     calEtaPt_.push_back(fs->make<TH1F>(cal.name(hist,"ktEtaPt",jdx), cal.name("calEtaPt",jdx), ((int)binsPt_.size()-1), &binsPt_[0]));
00157     resEtaPt_.push_back(fs->make<TH1F>(res.name(hist,"ktEtaPt",jdx), res.name("resEtaPt",jdx), ((int)binsPt_.size()-1), &binsPt_[0]));
00158     for(int idx=0; idx<((int)binsPt_.size()-1); ++idx)
00159       buffer.push_back( fs->make<TH1F>(fit.name(hist, "ktEtaPt",jdx,idx), fit.name("ktEtaPt",jdx,idx), bins_, min_, max_) );
00160     ktEtaPt_.push_back(buffer); 
00161   }
00162 }
00163 
00164 #endif