CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/RecoJets/JetAnalyzers/src/FFTJetPileupAnalyzer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetPileupAnalyzer
00004 // Class:      FFTJetPileupAnalyzer
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Thu Apr 21 15:52:11 CDT 2011
00016 // $Id: FFTJetPileupAnalyzer.cc,v 1.12 2011/07/18 17:40:54 igv Exp $
00017 //
00018 //
00019 
00020 #include <cassert>
00021 #include <sstream>
00022 #include <numeric>
00023 
00024 #include "TNtuple.h"
00025 
00026 // user include files
00027 #include "FWCore/Framework/interface/Frameworkfwd.h"
00028 #include "FWCore/Framework/interface/EDAnalyzer.h"
00029 #include "FWCore/Framework/interface/Event.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 #include "FWCore/ServiceRegistry/interface/Service.h"
00032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00033 
00034 #include "DataFormats/Common/interface/Handle.h"
00035 #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
00036 #include "DataFormats/JetReco/interface/FFTJetPileupSummary.h"
00037 #include "DataFormats/VertexReco/interface/Vertex.h"
00038 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00039 
00040 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
00041 
00042 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00043 
00044 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00045 
00046 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00047 
00048 //
00049 // class declaration
00050 //
00051 class FFTJetPileupAnalyzer : public edm::EDAnalyzer
00052 {
00053 public:
00054     explicit FFTJetPileupAnalyzer(const edm::ParameterSet&);
00055     ~FFTJetPileupAnalyzer();
00056 
00057 private:
00058     FFTJetPileupAnalyzer();
00059     FFTJetPileupAnalyzer(const FFTJetPileupAnalyzer&);
00060     FFTJetPileupAnalyzer& operator=(const FFTJetPileupAnalyzer&);
00061 
00062     // The following method should take all necessary info from
00063     // PileupSummaryInfo and fill out the ntuple
00064     void analyzePileup(const std::vector<PileupSummaryInfo>& pInfo);
00065 
00066     virtual void beginJob() ;
00067     virtual void analyze(const edm::Event&, const edm::EventSetup&);
00068     virtual void endJob() ;
00069 
00070     edm::InputTag histoLabel;
00071     edm::InputTag summaryLabel;
00072     edm::InputTag fastJetRhoLabel;
00073     edm::InputTag fastJetSigmaLabel;
00074     edm::InputTag gridLabel;
00075     edm::InputTag srcPVs;
00076     std::string pileupLabel;
00077     std::string ntupleName;
00078     std::string ntupleTitle;
00079     bool collectHistos;
00080     bool collectSummaries;
00081     bool collectFastJetRho;
00082     bool collectPileup;
00083     bool collectOOTPileup;
00084     bool collectGrids;
00085     bool collectGridDensity;
00086     bool collectVertexInfo;
00087     bool verbosePileupInfo;
00088 
00089     double vertexNdofCut;
00090     double crazyEnergyCut;
00091 
00092     std::vector<float> ntupleData;
00093     TNtuple* nt;
00094     int totalNpu;
00095     int totalNPV;
00096     unsigned long counter;
00097 };
00098 
00099 //
00100 // constructors and destructor
00101 //
00102 FFTJetPileupAnalyzer::FFTJetPileupAnalyzer(const edm::ParameterSet& ps)
00103     : init_param(edm::InputTag, histoLabel),
00104       init_param(edm::InputTag, summaryLabel),
00105       init_param(edm::InputTag, fastJetRhoLabel),
00106       init_param(edm::InputTag, fastJetSigmaLabel),
00107       init_param(edm::InputTag, gridLabel),
00108       init_param(edm::InputTag, srcPVs),
00109       init_param(std::string, pileupLabel),
00110       init_param(std::string, ntupleName),
00111       init_param(std::string, ntupleTitle),
00112       init_param(bool, collectHistos),
00113       init_param(bool, collectSummaries),
00114       init_param(bool, collectFastJetRho),
00115       init_param(bool, collectPileup),
00116       init_param(bool, collectOOTPileup),
00117       init_param(bool, collectGrids),
00118       init_param(bool, collectGridDensity),
00119       init_param(bool, collectVertexInfo),
00120       init_param(bool, verbosePileupInfo),
00121       init_param(double, vertexNdofCut),
00122       init_param(double, crazyEnergyCut),
00123       nt(0),
00124       totalNpu(-1),
00125       totalNPV(-1),
00126       counter(0)
00127 {
00128 }
00129 
00130 
00131 FFTJetPileupAnalyzer::~FFTJetPileupAnalyzer()
00132 {
00133 }
00134 
00135 
00136 //
00137 // member functions
00138 //
00139 void FFTJetPileupAnalyzer::analyzePileup(
00140     const std::vector<PileupSummaryInfo>& info)
00141 {
00142     const unsigned nBx = info.size();
00143     if (collectPileup)
00144         ntupleData.push_back(static_cast<float>(nBx));
00145 
00146     double sumpt_Lo = 0.0, sumpt_Hi = 0.0;
00147     totalNpu = 0;
00148 
00149     int npu_by_Bx[3] = {0,};
00150     double sumpt_Lo_by_Bx[3] = {0.0,}, sumpt_Hi_by_Bx[3] = {0.0,};
00151 
00152     if (verbosePileupInfo)
00153         std::cout << "\n**** Pileup info begin" << std::endl;
00154 
00155     bool isCrazy = false;
00156     for (unsigned ibx = 0; ibx < nBx; ++ibx)
00157     {
00158         const PileupSummaryInfo& puInfo(info[ibx]);
00159 
00160         const int bx = puInfo.getBunchCrossing();
00161         const int npu = puInfo.getPU_NumInteractions();
00162         const std::vector<float>& lopt(puInfo.getPU_sumpT_lowpT());
00163         const std::vector<float>& hipt(puInfo.getPU_sumpT_highpT());
00164         const double losum = std::accumulate(lopt.begin(), lopt.end(), 0.0);
00165         const double hisum = std::accumulate(hipt.begin(), hipt.end(), 0.0);
00166 
00167         if (losum >= crazyEnergyCut)
00168             isCrazy = true;
00169         if (hisum >= crazyEnergyCut)
00170             isCrazy = true;
00171 
00172         totalNpu += npu;
00173         sumpt_Lo += losum;
00174         sumpt_Hi += hisum;
00175 
00176         const unsigned idx = bx < 0 ? 0U : (bx == 0 ? 1U : 2U);
00177         npu_by_Bx[idx] += npu;
00178         sumpt_Lo_by_Bx[idx] += losum;
00179         sumpt_Hi_by_Bx[idx] += hisum;
00180 
00181         if (verbosePileupInfo)
00182             std::cout << "ibx " << ibx << " bx " << bx
00183                       << " npu " << npu << " losum " << losum
00184                       << " hisum " << hisum
00185                       << std::endl;
00186     }
00187 
00188     if (verbosePileupInfo)
00189         std::cout << "**** Pileup info end\n" << std::endl;
00190 
00191     if (isCrazy)
00192     {
00193         totalNpu = -1;
00194         sumpt_Lo = 0.0;
00195         sumpt_Hi = 0.0;
00196         for (unsigned ibx = 0; ibx < 3; ++ibx)
00197         {
00198             npu_by_Bx[ibx] = -1;
00199             sumpt_Lo_by_Bx[ibx] = 0.0;
00200             sumpt_Hi_by_Bx[ibx] = 0.0;
00201         }
00202     }
00203 
00204     if (collectPileup)
00205     {
00206         ntupleData.push_back(totalNpu);
00207         ntupleData.push_back(sumpt_Lo);
00208         ntupleData.push_back(sumpt_Hi);
00209     }
00210 
00211     if (collectOOTPileup)
00212         for (unsigned ibx = 0; ibx < 3; ++ibx)
00213         {
00214             ntupleData.push_back(npu_by_Bx[ibx]);
00215             ntupleData.push_back(sumpt_Lo_by_Bx[ibx]);
00216             ntupleData.push_back(sumpt_Hi_by_Bx[ibx]);
00217         }
00218 }
00219 
00220 
00221 // ------------ method called once each job just before starting event loop
00222 void FFTJetPileupAnalyzer::beginJob()
00223 {
00224     // Come up with the list of variables
00225     std::string vars = "cnt:run:event";
00226     if (collectPileup)
00227         vars += ":nbx:npu:sumptLowCut:sumptHiCut";
00228     if (collectOOTPileup)
00229     {
00230         vars += ":npu_negbx:sumptLowCut_negbx:sumptHiCut_negbx";
00231         vars += ":npu_0bx:sumptLowCut_0bx:sumptHiCut_0bx";
00232         vars += ":npu_posbx:sumptLowCut_posbx:sumptHiCut_posbx";
00233     }
00234     if (collectSummaries)
00235         vars += ":estimate:pileup:uncert:uncertCode";
00236     if (collectFastJetRho)
00237         vars += ":fjrho:fjsigma";
00238     if (collectGridDensity)
00239         vars += ":gridEtDensity:gridEtDensityMixed";
00240     if (collectVertexInfo)
00241         vars += ":nPV";
00242 
00243     // Book the ntuple
00244     edm::Service<TFileService> fs;
00245     nt = fs->make<TNtuple>(ntupleName.c_str(), ntupleTitle.c_str(),
00246                            vars.c_str());
00247     ntupleData.reserve(nt->GetNvar());
00248 }
00249 
00250 
00251 // ------------ method called to for each event  ------------
00252 void FFTJetPileupAnalyzer::analyze(const edm::Event& iEvent,
00253                                    const edm::EventSetup& iSetup)
00254 {
00255     ntupleData.clear();
00256     ntupleData.push_back(counter);
00257     totalNpu = -1;
00258     totalNPV = -1;
00259 
00260     const long runnumber = iEvent.id().run();
00261     const long eventnumber = iEvent.id().event();
00262     ntupleData.push_back(runnumber);
00263     ntupleData.push_back(eventnumber);
00264 
00265     // Get pileup information from the pile-up information module
00266     if (collectPileup || collectOOTPileup)
00267     {
00268         edm::Handle<std::vector<PileupSummaryInfo> > puInfo;
00269         if (iEvent.getByLabel(pileupLabel, puInfo))
00270             analyzePileup(*puInfo);
00271         else
00272         {
00273             if (collectPileup)
00274             {
00275                 ntupleData.push_back(-1);
00276                 ntupleData.push_back(-1);
00277                 ntupleData.push_back(0.f);
00278                 ntupleData.push_back(0.f);
00279             }
00280             if (collectOOTPileup)
00281                 for (unsigned ibx = 0; ibx < 3; ++ibx)
00282                 {
00283                     ntupleData.push_back(-1);
00284                     ntupleData.push_back(0.f);
00285                     ntupleData.push_back(0.f);
00286                 }
00287         }
00288     }
00289 
00290     if (collectHistos)
00291     {
00292         edm::Handle<TH2D> input;
00293         iEvent.getByLabel(histoLabel, input);
00294 
00295         edm::Service<TFileService> fs;
00296         TH2D* copy = new TH2D(*input);
00297 
00298         std::ostringstream os;
00299         os << copy->GetName() << '_' << counter << '_'
00300            << totalNpu << '_' << runnumber << '_' << eventnumber;
00301         const std::string& newname(os.str());
00302         copy->SetNameTitle(newname.c_str(), newname.c_str());
00303 
00304         copy->SetDirectory(fs->getBareDirectory());
00305     }
00306 
00307     if (collectSummaries)
00308     {
00309         edm::Handle<reco::FFTJetPileupSummary> summary;
00310         iEvent.getByLabel(summaryLabel, summary);
00311 
00312         ntupleData.push_back(summary->uncalibratedQuantile());
00313         ntupleData.push_back(summary->pileupRho());
00314         ntupleData.push_back(summary->pileupRhoUncertainty());
00315         ntupleData.push_back(summary->uncertaintyCode());
00316     }
00317 
00318     if (collectFastJetRho)
00319     {
00320         edm::Handle<double> fjrho, fjsigma;
00321         iEvent.getByLabel(fastJetRhoLabel, fjrho);
00322         iEvent.getByLabel(fastJetSigmaLabel, fjsigma);
00323 
00324         ntupleData.push_back(*fjrho);
00325         ntupleData.push_back(*fjsigma);
00326     }
00327 
00328     if (collectGrids)
00329     {
00330         edm::Handle<reco::DiscretizedEnergyFlow> input;
00331         iEvent.getByLabel(gridLabel, input);
00332 
00333         // Make sure the input grid is reasonable
00334         const double* data = input->data();
00335         assert(data);
00336         assert(input->phiBin0Edge() == 0.0);
00337         const unsigned nEta = input->nEtaBins();
00338         const unsigned nPhi = input->nPhiBins();
00339 
00340         // Generate a name for the output histogram
00341         std::ostringstream os;
00342         os << "FFTJetGrid_" << counter << '_'
00343            << totalNpu << '_' << runnumber << '_' << eventnumber;
00344         const std::string& newname(os.str());
00345 
00346         // Make a histogram and copy the grid data into it
00347         edm::Service<TFileService> fs;
00348         TH2F* h = fs->make<TH2F>(newname.c_str(), newname.c_str(),
00349                                  nEta, input->etaMin(), input->etaMax(),
00350                                  nPhi, 0.0, 2.0*M_PI);
00351         h->GetXaxis()->SetTitle("Eta");
00352         h->GetYaxis()->SetTitle("Phi");
00353         h->GetZaxis()->SetTitle("Transverse Energy");
00354 
00355         for (unsigned ieta=0; ieta<nEta; ++ieta)
00356             for (unsigned iphi=0; iphi<nPhi; ++iphi)
00357                 h->SetBinContent(ieta+1U, iphi+1U, data[ieta*nPhi + iphi]);
00358     }
00359 
00360     if (collectGridDensity)
00361     {
00362         edm::Handle<std::pair<double,double> > etSum;
00363         iEvent.getByLabel(histoLabel, etSum);
00364 
00365         ntupleData.push_back(etSum->first);
00366         ntupleData.push_back(etSum->second);
00367     }
00368 
00369     if (collectVertexInfo)
00370     {
00371         edm::Handle<reco::VertexCollection> pvCollection;
00372         iEvent.getByLabel(srcPVs, pvCollection);
00373         totalNPV = 0;
00374         if (!pvCollection->empty())
00375             for (reco::VertexCollection::const_iterator pv = pvCollection->begin();
00376                  pv != pvCollection->end(); ++pv)
00377             {
00378                 const double ndof = pv->ndof();
00379                 if (!pv->isFake() && ndof > vertexNdofCut)
00380                     ++totalNPV;
00381             }
00382         ntupleData.push_back(totalNPV);
00383     }
00384 
00385     assert(ntupleData.size() == static_cast<unsigned>(nt->GetNvar()));
00386     nt->Fill(&ntupleData[0]);
00387 
00388     ++counter;
00389 }
00390 
00391 
00392 // ------------ method called once each job just after ending the event loop
00393 void FFTJetPileupAnalyzer::endJob()
00394 {
00395 }
00396 
00397 
00398 //define this as a plug-in
00399 DEFINE_FWK_MODULE(FFTJetPileupAnalyzer);