#include <RecoJets/FFTJetProducers/plugins/FFTJetPileupEstimator.cc>
Public Member Functions | |
FFTJetPileupEstimator (const edm::ParameterSet &) | |
~FFTJetPileupEstimator () | |
Protected Member Functions | |
void | beginJob () |
void | endJob () |
void | produce (edm::Event &, const edm::EventSetup &) |
Private Member Functions | |
template<class Ptr > | |
void | checkConfig (const Ptr &ptr, const char *message) |
FFTJetPileupEstimator (const FFTJetPileupEstimator &) | |
FFTJetPileupEstimator () | |
FFTJetPileupEstimator & | operator= (const FFTJetPileupEstimator &) |
Private Attributes | |
std::auto_ptr < fftjet::Functor1< double, double > > | calibrationCurve |
double | cdfvalue |
unsigned | filterNumber |
edm::InputTag | inputLabel |
std::string | outputLabel |
double | ptToDensityFactor |
std::auto_ptr < fftjet::Functor1< double, double > > | uncertaintyCurve |
std::vector< double > | uncertaintyZones |
Description: estimates the actual pileup
Implementation: [Notes on implementation]
Definition at line 43 of file FFTJetPileupEstimator.cc.
FFTJetPileupEstimator::FFTJetPileupEstimator | ( | const edm::ParameterSet & | ps | ) | [explicit] |
Definition at line 80 of file FFTJetPileupEstimator.cc.
References calibrationCurve, checkConfig(), fftjetcms::fftjet_Function_parser(), edm::ParameterSet::getParameter(), outputLabel, and uncertaintyCurve.
: init_param(edm::InputTag, inputLabel), init_param(std::string, outputLabel), init_param(double, cdfvalue), init_param(double, ptToDensityFactor), init_param(unsigned, filterNumber), init_param(std::vector<double>, uncertaintyZones) { calibrationCurve = fftjet_Function_parser( ps.getParameter<edm::ParameterSet>("calibrationCurve")); checkConfig(calibrationCurve, "bad calibration curve definition"); uncertaintyCurve = fftjet_Function_parser( ps.getParameter<edm::ParameterSet>("uncertaintyCurve")); checkConfig(uncertaintyCurve, "bad uncertainty curve definition"); produces<reco::FFTJetPileupSummary>(outputLabel); }
FFTJetPileupEstimator::~FFTJetPileupEstimator | ( | ) |
Definition at line 100 of file FFTJetPileupEstimator.cc.
{ }
FFTJetPileupEstimator::FFTJetPileupEstimator | ( | ) | [private] |
FFTJetPileupEstimator::FFTJetPileupEstimator | ( | const FFTJetPileupEstimator & | ) | [private] |
void FFTJetPileupEstimator::beginJob | ( | void | ) | [protected, virtual] |
void FFTJetPileupEstimator::checkConfig | ( | const Ptr & | ptr, |
const char * | message | ||
) | [inline, private] |
Definition at line 61 of file FFTJetPileupEstimator.cc.
References NULL.
Referenced by FFTJetPileupEstimator().
{ if (ptr.get() == NULL) throw cms::Exception("FFTJetBadConfig") << message << std::endl; }
void FFTJetPileupEstimator::endJob | ( | void | ) | [protected, virtual] |
FFTJetPileupEstimator& FFTJetPileupEstimator::operator= | ( | const FFTJetPileupEstimator & | ) | [private] |
void FFTJetPileupEstimator::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [protected, virtual] |
Implements edm::EDProducer.
Definition at line 109 of file FFTJetPileupEstimator.cc.
References cdfvalue, Exception, filterNumber, edm::Event::getByLabel(), h, i, collect_tpl::input, inputLabel, outputLabel, ptToDensityFactor, edm::Event::put(), summarizeEdmComparisonLogfiles::summary, and uncertaintyZones.
{ edm::Handle<TH2D> input; iEvent.getByLabel(inputLabel, input); const TH2D& h(*input); const unsigned nScales = h.GetXaxis()->GetNbins(); const unsigned nCdfvalues = h.GetYaxis()->GetNbins(); const unsigned fixedCdfvalueBin = static_cast<unsigned>( std::floor(cdfvalue*nCdfvalues)); if (fixedCdfvalueBin >= nCdfvalues) { throw cms::Exception("FFTJetBadConfig") << "Bad cdf value" << std::endl; } if (filterNumber >= nScales) { throw cms::Exception("FFTJetBadConfig") << "Bad filter number" << std::endl; } // Simple fixed-point pile-up estimate const double curve = h.GetBinContent(filterNumber+1U, fixedCdfvalueBin+1U); const double pileupRho = ptToDensityFactor*(*calibrationCurve)(curve); const double rhoUncert = ptToDensityFactor*(*uncertaintyCurve)(curve); // Determine the uncertainty zone of the estimate. The "curve" // has to be above or equal to uncertaintyZones[i] but below // uncertaintyZones[i + 1] (the second condition is also satisfied // by i == uncertaintyZones.size() - 1). Of course, it is assumed // that the vector of zones is configured appropriately -- the zone // boundaries must be presented in the increasing order. int uncertaintyCode = -1; if (!uncertaintyZones.empty()) { const unsigned nZones = uncertaintyZones.size(); for (unsigned i = 0; i < nZones; ++i) if (curve >= uncertaintyZones[i]) { if (i == nZones - 1U) { uncertaintyCode = i; break; } else if (curve < uncertaintyZones[i + 1]) { uncertaintyCode = i; break; } } } std::auto_ptr<reco::FFTJetPileupSummary> summary( new reco::FFTJetPileupSummary(curve, pileupRho, rhoUncert, uncertaintyCode)); iEvent.put(summary, outputLabel); }
std::auto_ptr<fftjet::Functor1<double,double> > FFTJetPileupEstimator::calibrationCurve [private] |
Definition at line 73 of file FFTJetPileupEstimator.cc.
Referenced by FFTJetPileupEstimator().
double FFTJetPileupEstimator::cdfvalue [private] |
Definition at line 69 of file FFTJetPileupEstimator.cc.
Referenced by produce().
unsigned FFTJetPileupEstimator::filterNumber [private] |
Definition at line 71 of file FFTJetPileupEstimator.cc.
Referenced by produce().
Definition at line 67 of file FFTJetPileupEstimator.cc.
Referenced by produce().
std::string FFTJetPileupEstimator::outputLabel [private] |
Definition at line 68 of file FFTJetPileupEstimator.cc.
Referenced by FFTJetPileupEstimator(), and produce().
double FFTJetPileupEstimator::ptToDensityFactor [private] |
Definition at line 70 of file FFTJetPileupEstimator.cc.
Referenced by produce().
std::auto_ptr<fftjet::Functor1<double,double> > FFTJetPileupEstimator::uncertaintyCurve [private] |
Definition at line 74 of file FFTJetPileupEstimator.cc.
Referenced by FFTJetPileupEstimator().
std::vector<double> FFTJetPileupEstimator::uncertaintyZones [private] |
Definition at line 72 of file FFTJetPileupEstimator.cc.
Referenced by produce().