00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020 #include <cmath>
00021
00022
00023 #include "FWCore/Framework/interface/EDProducer.h"
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/Framework/interface/MakerMacros.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027
00028
00029 #include "DataFormats/Common/interface/View.h"
00030 #include "DataFormats/Common/interface/Handle.h"
00031
00032 #include "DataFormats/JetReco/interface/FFTJetPileupSummary.h"
00033 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00034
00035 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00036
00037
00038 #include "JetMETCorrections/FFTJetObjects/interface/FFTJetLookupTableSequenceLoader.h"
00039
00040 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00041
00042 using namespace fftjetcms;
00043
00044
00045
00046
00047 class FFTJetPileupEstimator : public edm::EDProducer
00048 {
00049 public:
00050 explicit FFTJetPileupEstimator(const edm::ParameterSet&);
00051 ~FFTJetPileupEstimator();
00052
00053 protected:
00054
00055 void beginJob();
00056 void produce(edm::Event&, const edm::EventSetup&);
00057 void endJob();
00058
00059 private:
00060 FFTJetPileupEstimator();
00061 FFTJetPileupEstimator(const FFTJetPileupEstimator&);
00062 FFTJetPileupEstimator& operator=(const FFTJetPileupEstimator&);
00063
00064 std::auto_ptr<reco::FFTJetPileupSummary> calibrateFromConfig(
00065 double uncalibrated) const;
00066
00067 std::auto_ptr<reco::FFTJetPileupSummary> calibrateFromDB(
00068 double uncalibrated, const edm::EventSetup& iSetup) const;
00069
00070 template<class Ptr>
00071 inline void checkConfig(const Ptr& ptr, const char* message)
00072 {
00073 if (ptr.get() == NULL)
00074 throw cms::Exception("FFTJetBadConfig") << message << std::endl;
00075 }
00076
00077 edm::InputTag inputLabel;
00078 std::string outputLabel;
00079 double cdfvalue;
00080 double ptToDensityFactor;
00081 unsigned filterNumber;
00082 std::vector<double> uncertaintyZones;
00083 std::auto_ptr<fftjet::Functor1<double,double> > calibrationCurve;
00084 std::auto_ptr<fftjet::Functor1<double,double> > uncertaintyCurve;
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099 std::string calibTableRecord;
00100 std::string calibTableCategory;
00101 std::string uncertaintyZonesName;
00102 std::string calibrationCurveName;
00103 std::string uncertaintyCurveName;
00104 bool loadCalibFromDB;
00105 };
00106
00107
00108
00109
00110 FFTJetPileupEstimator::FFTJetPileupEstimator(const edm::ParameterSet& ps)
00111 : init_param(edm::InputTag, inputLabel),
00112 init_param(std::string, outputLabel),
00113 init_param(double, cdfvalue),
00114 init_param(double, ptToDensityFactor),
00115 init_param(unsigned, filterNumber),
00116 init_param(std::vector<double>, uncertaintyZones),
00117 init_param(std::string, calibTableRecord),
00118 init_param(std::string, calibTableCategory),
00119 init_param(std::string, uncertaintyZonesName),
00120 init_param(std::string, calibrationCurveName),
00121 init_param(std::string, uncertaintyCurveName),
00122 init_param(bool, loadCalibFromDB)
00123 {
00124 calibrationCurve = fftjet_Function_parser(
00125 ps.getParameter<edm::ParameterSet>("calibrationCurve"));
00126 checkConfig(calibrationCurve, "bad calibration curve definition");
00127
00128 uncertaintyCurve = fftjet_Function_parser(
00129 ps.getParameter<edm::ParameterSet>("uncertaintyCurve"));
00130 checkConfig(uncertaintyCurve, "bad uncertainty curve definition");
00131
00132 produces<reco::FFTJetPileupSummary>(outputLabel);
00133 }
00134
00135
00136 FFTJetPileupEstimator::~FFTJetPileupEstimator()
00137 {
00138 }
00139
00140
00141
00142
00143
00144
00145 void FFTJetPileupEstimator::produce(edm::Event& iEvent,
00146 const edm::EventSetup& iSetup)
00147 {
00148 edm::Handle<reco::DiscretizedEnergyFlow> input;
00149 iEvent.getByLabel(inputLabel, input);
00150
00151 const reco::DiscretizedEnergyFlow& h(*input);
00152 const unsigned nScales = h.nEtaBins();
00153 const unsigned nCdfvalues = h.nPhiBins();
00154
00155 const unsigned fixedCdfvalueBin = static_cast<unsigned>(
00156 std::floor(cdfvalue*nCdfvalues));
00157 if (fixedCdfvalueBin >= nCdfvalues)
00158 {
00159 throw cms::Exception("FFTJetBadConfig")
00160 << "Bad cdf value" << std::endl;
00161 }
00162 if (filterNumber >= nScales)
00163 {
00164 throw cms::Exception("FFTJetBadConfig")
00165 << "Bad filter number" << std::endl;
00166 }
00167
00168
00169 const double curve = h.data()[filterNumber*nCdfvalues + fixedCdfvalueBin];
00170
00171 std::auto_ptr<reco::FFTJetPileupSummary> summary;
00172 if (loadCalibFromDB)
00173 summary = calibrateFromDB(curve, iSetup);
00174 else
00175 summary = calibrateFromConfig(curve);
00176 iEvent.put(summary, outputLabel);
00177 }
00178
00179
00180 void FFTJetPileupEstimator::beginJob()
00181 {
00182 }
00183
00184
00185 void FFTJetPileupEstimator::endJob()
00186 {
00187 }
00188
00189
00190 std::auto_ptr<reco::FFTJetPileupSummary>
00191 FFTJetPileupEstimator::calibrateFromConfig(const double curve) const
00192 {
00193 const double pileupRho = ptToDensityFactor*(*calibrationCurve)(curve);
00194 const double rhoUncert = ptToDensityFactor*(*uncertaintyCurve)(curve);
00195
00196
00197
00198
00199
00200
00201
00202 int uncertaintyCode = -1;
00203 if (!uncertaintyZones.empty())
00204 {
00205 const unsigned nZones = uncertaintyZones.size();
00206 for (unsigned i = 0; i < nZones; ++i)
00207 if (curve >= uncertaintyZones[i])
00208 {
00209 if (i == nZones - 1U)
00210 {
00211 uncertaintyCode = i;
00212 break;
00213 }
00214 else if (curve < uncertaintyZones[i + 1])
00215 {
00216 uncertaintyCode = i;
00217 break;
00218 }
00219 }
00220 }
00221
00222 return std::auto_ptr<reco::FFTJetPileupSummary>(
00223 new reco::FFTJetPileupSummary(curve, pileupRho,
00224 rhoUncert, uncertaintyCode));
00225 }
00226
00227
00228 std::auto_ptr<reco::FFTJetPileupSummary>
00229 FFTJetPileupEstimator::calibrateFromDB(
00230 const double curve, const edm::EventSetup& iSetup) const
00231 {
00232 edm::ESHandle<FFTJetLookupTableSequence> h;
00233 StaticFFTJetLookupTableSequenceLoader::instance().load(
00234 iSetup, calibTableRecord, h);
00235 boost::shared_ptr<npstat::StorableMultivariateFunctor> uz =
00236 (*h)[calibTableCategory][uncertaintyZonesName];
00237 boost::shared_ptr<npstat::StorableMultivariateFunctor> cc =
00238 (*h)[calibTableCategory][calibrationCurveName];
00239 boost::shared_ptr<npstat::StorableMultivariateFunctor> uc =
00240 (*h)[calibTableCategory][uncertaintyCurveName];
00241
00242 const double pileupRho = ptToDensityFactor*(*cc)(&curve, 1U);
00243 const double rhoUncert = ptToDensityFactor*(*uc)(&curve, 1U);
00244 const int uncertaintyCode = round((*uz)(&curve, 1U));
00245
00246 return std::auto_ptr<reco::FFTJetPileupSummary>(
00247 new reco::FFTJetPileupSummary(curve, pileupRho,
00248 rhoUncert, uncertaintyCode));
00249 }
00250
00251
00252
00253 DEFINE_FWK_MODULE(FFTJetPileupEstimator);