00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020 #include <cmath>
00021 #include <fstream>
00022
00023
00024 #include "fftjet/FrequencyKernelConvolver.hh"
00025 #include "fftjet/DiscreteGauss2d.hh"
00026 #include "fftjet/EquidistantSequence.hh"
00027 #include "fftjet/interpolate.hh"
00028
00029
00030 #include "FWCore/Framework/interface/Frameworkfwd.h"
00031 #include "FWCore/Framework/interface/EDProducer.h"
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034
00035 #include "DataFormats/Common/interface/View.h"
00036 #include "DataFormats/Common/interface/Handle.h"
00037
00038 #include "DataFormats/JetReco/interface/DiscretizedEnergyFlow.h"
00039
00040 #include "RecoJets/FFTJetAlgorithms/interface/gridConverters.h"
00041
00042
00043 #include "RecoJets/FFTJetProducers/interface/FFTJetParameterParser.h"
00044
00045
00046 #include "RecoJets/FFTJetProducers/interface/FFTJetInterface.h"
00047
00048
00049 #include "JetMETCorrections/FFTJetObjects/interface/FFTJetLookupTableSequenceLoader.h"
00050
00051
00052 using namespace fftjetcms;
00053
00054
00055
00056
00057 class FFTJetPileupProcessor : public edm::EDProducer, public FFTJetInterface
00058 {
00059 public:
00060 explicit FFTJetPileupProcessor(const edm::ParameterSet&);
00061 ~FFTJetPileupProcessor();
00062
00063 protected:
00064
00065 void beginJob() ;
00066 void produce(edm::Event&, const edm::EventSetup&);
00067 void endJob() ;
00068
00069 private:
00070 FFTJetPileupProcessor();
00071 FFTJetPileupProcessor(const FFTJetPileupProcessor&);
00072 FFTJetPileupProcessor& operator=(const FFTJetPileupProcessor&);
00073
00074 void buildKernelConvolver(const edm::ParameterSet&);
00075 void mixExtraGrid();
00076 void loadFlatteningFactors(const edm::EventSetup& iSetup);
00077
00078
00079 std::auto_ptr<MyFFTEngine> engine;
00080
00081
00082 std::auto_ptr<fftjet::AbsFrequencyKernel> kernel2d;
00083
00084
00085 std::auto_ptr<fftjet::AbsConvolverBase<Real> > convolver;
00086
00087
00088 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > convolvedFlow;
00089
00090
00091 std::auto_ptr<fftjet::EquidistantInLogSpace> filterScales;
00092
00093
00094
00095 std::vector<double> etaFlatteningFactors;
00096
00097
00098 unsigned nPercentiles;
00099
00100
00101 unsigned convolverMinBin;
00102 unsigned convolverMaxBin;
00103
00104
00105 double pileupEtaPhiArea;
00106
00107
00108 std::vector<std::string> externalGridFiles;
00109 std::ifstream gridStream;
00110 double externalGridMaxEnergy;
00111 unsigned currentFileNum;
00112
00113
00114 std::vector<double> percentileData;
00115
00116
00117
00118
00119
00120 std::string flatteningTableRecord;
00121 std::string flatteningTableName;
00122 std::string flatteningTableCategory;
00123 bool loadFlatteningFactorsFromDB;
00124 };
00125
00126
00127
00128
00129 FFTJetPileupProcessor::FFTJetPileupProcessor(const edm::ParameterSet& ps)
00130 : FFTJetInterface(ps),
00131 etaFlatteningFactors(
00132 ps.getParameter<std::vector<double> >("etaFlatteningFactors")),
00133 nPercentiles(ps.getParameter<unsigned>("nPercentiles")),
00134 convolverMinBin(ps.getParameter<unsigned>("convolverMinBin")),
00135 convolverMaxBin(ps.getParameter<unsigned>("convolverMaxBin")),
00136 pileupEtaPhiArea(ps.getParameter<double>("pileupEtaPhiArea")),
00137 externalGridFiles(ps.getParameter<std::vector<std::string> >("externalGridFiles")),
00138 externalGridMaxEnergy(ps.getParameter<double>("externalGridMaxEnergy")),
00139 currentFileNum(externalGridFiles.size() + 1U),
00140 flatteningTableRecord(ps.getParameter<std::string>("flatteningTableRecord")),
00141 flatteningTableName(ps.getParameter<std::string>("flatteningTableName")),
00142 flatteningTableCategory(ps.getParameter<std::string>("flatteningTableCategory")),
00143 loadFlatteningFactorsFromDB(ps.getParameter<bool>("loadFlatteningFactorsFromDB"))
00144 {
00145
00146 energyFlow = fftjet_Grid2d_parser(
00147 ps.getParameter<edm::ParameterSet>("GridConfiguration"));
00148 checkConfig(energyFlow, "invalid discretization grid");
00149
00150
00151 convolvedFlow = std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >(
00152 new fftjet::Grid2d<fftjetcms::Real>(*energyFlow));
00153
00154
00155 if (!etaFlatteningFactors.empty())
00156 {
00157 if (etaFlatteningFactors.size() != convolvedFlow->nEta())
00158 throw cms::Exception("FFTJetBadConfig")
00159 << "ERROR in FFTJetPileupProcessor constructor:"
00160 " number of elements in the \"etaFlatteningFactors\""
00161 " vector is inconsistent with the discretization grid binning"
00162 << std::endl;
00163 }
00164
00165
00166
00167 buildKernelConvolver(ps);
00168
00169
00170 const unsigned nScales = ps.getParameter<unsigned>("nScales");
00171 const double minScale = ps.getParameter<double>("minScale");
00172 const double maxScale = ps.getParameter<double>("maxScale");
00173 if (minScale <= 0.0 || maxScale < minScale || nScales == 0U)
00174 throw cms::Exception("FFTJetBadConfig")
00175 << "invalid filter scales" << std::endl;
00176
00177 filterScales = std::auto_ptr<fftjet::EquidistantInLogSpace>(
00178 new fftjet::EquidistantInLogSpace(minScale, maxScale, nScales));
00179
00180 percentileData.resize(nScales*nPercentiles);
00181
00182 produces<reco::DiscretizedEnergyFlow>(outputLabel);
00183 produces<std::pair<double,double> >(outputLabel);
00184 }
00185
00186
00187 void FFTJetPileupProcessor::buildKernelConvolver(const edm::ParameterSet& ps)
00188 {
00189
00190 double kernelEtaScale = ps.getParameter<double>("kernelEtaScale");
00191 const double kernelPhiScale = ps.getParameter<double>("kernelPhiScale");
00192 if (kernelEtaScale <= 0.0 || kernelPhiScale <= 0.0)
00193 throw cms::Exception("FFTJetBadConfig")
00194 << "invalid kernel scales" << std::endl;
00195
00196 if (convolverMinBin || convolverMaxBin)
00197 if (convolverMinBin >= convolverMaxBin ||
00198 convolverMaxBin > energyFlow->nEta())
00199 throw cms::Exception("FFTJetBadConfig")
00200 << "invalid convolver bin range" << std::endl;
00201
00202
00203
00204 kernelEtaScale *= (2.0*M_PI/(energyFlow->etaMax() - energyFlow->etaMin()));
00205
00206
00207 engine = std::auto_ptr<MyFFTEngine>(
00208 new MyFFTEngine(energyFlow->nEta(), energyFlow->nPhi()));
00209
00210
00211 kernel2d = std::auto_ptr<fftjet::AbsFrequencyKernel>(
00212 new fftjet::DiscreteGauss2d(
00213 kernelEtaScale, kernelPhiScale,
00214 energyFlow->nEta(), energyFlow->nPhi()));
00215
00216
00217 convolver = std::auto_ptr<fftjet::AbsConvolverBase<Real> >(
00218 new fftjet::FrequencyKernelConvolver<Real,Complex>(
00219 engine.get(), kernel2d.get(),
00220 convolverMinBin, convolverMaxBin));
00221 }
00222
00223
00224 FFTJetPileupProcessor::~FFTJetPileupProcessor()
00225 {
00226 }
00227
00228
00229
00230 void FFTJetPileupProcessor::produce(
00231 edm::Event& iEvent, const edm::EventSetup& iSetup)
00232 {
00233 loadInputCollection(iEvent);
00234 discretizeEnergyFlow();
00235
00236
00237
00238 const fftjet::Grid2d<Real>& g(*energyFlow);
00239 const double densityBeforeMixing = g.sum()/pileupEtaPhiArea;
00240
00241
00242 double densityAfterMixing = -1.0;
00243 if (!externalGridFiles.empty())
00244 {
00245 mixExtraGrid();
00246 densityAfterMixing = g.sum()/pileupEtaPhiArea;
00247 }
00248
00249
00250 const unsigned nScales = filterScales->size();
00251 const double* scales = &(*filterScales)[0];
00252 Real* convData = const_cast<Real*>(convolvedFlow->data());
00253 Real* sortData = convData + convolverMinBin*convolvedFlow->nPhi();
00254 const unsigned dataLen = convolverMaxBin ?
00255 (convolverMaxBin - convolverMinBin)*convolvedFlow->nPhi() :
00256 convolvedFlow->nEta()*convolvedFlow->nPhi();
00257
00258
00259 if (loadFlatteningFactorsFromDB)
00260 loadFlatteningFactors(iSetup);
00261
00262
00263 convolver->setEventData(g.data(), g.nEta(), g.nPhi());
00264 for (unsigned iscale=0; iscale<nScales; ++iscale)
00265 {
00266
00267 convolver->convolveWithKernel(
00268 scales[iscale], convData,
00269 convolvedFlow->nEta(), convolvedFlow->nPhi());
00270
00271
00272 if (!etaFlatteningFactors.empty())
00273 convolvedFlow->scaleData(&etaFlatteningFactors[0],
00274 etaFlatteningFactors.size());
00275
00276
00277 std::sort(sortData, sortData+dataLen);
00278
00279
00280 for (unsigned iper=0; iper<nPercentiles; ++iper)
00281 {
00282
00283
00284 const double q = (iper + 0.5)/nPercentiles;
00285 const double dindex = q*(dataLen-1U);
00286 const unsigned ilow = static_cast<unsigned>(std::floor(dindex));
00287 const double percentile = fftjet::lin_interpolate_1d(
00288 ilow, ilow+1U, sortData[ilow], sortData[ilow+1U], dindex);
00289
00290
00291 percentileData[iscale*nPercentiles + iper] = percentile;
00292 }
00293 }
00294
00295
00296
00297 std::auto_ptr<reco::DiscretizedEnergyFlow> pTable(
00298 new reco::DiscretizedEnergyFlow(
00299 &percentileData[0], "FFTJetPileupProcessor",
00300 -0.5, nScales-0.5, 0.0, nScales, nPercentiles));
00301 iEvent.put(pTable, outputLabel);
00302
00303 std::auto_ptr<std::pair<double,double> > etSum(
00304 new std::pair<double,double>(densityBeforeMixing, densityAfterMixing));
00305 iEvent.put(etSum, outputLabel);
00306 }
00307
00308
00309 void FFTJetPileupProcessor::mixExtraGrid()
00310 {
00311 const unsigned nFiles = externalGridFiles.size();
00312 if (currentFileNum > nFiles)
00313 {
00314
00315 currentFileNum = 0;
00316 gridStream.open(externalGridFiles[currentFileNum].c_str(),
00317 std::ios_base::in | std::ios_base::binary);
00318 if (!gridStream.is_open())
00319 throw cms::Exception("FFTJetBadConfig")
00320 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00321 " failed to open external grid file "
00322 << externalGridFiles[currentFileNum] << std::endl;
00323 }
00324
00325 const fftjet::Grid2d<float>* g = 0;
00326 const unsigned maxFail = 100U;
00327 unsigned nEnergyRejected = 0;
00328
00329 while(!g)
00330 {
00331 g = fftjet::Grid2d<float>::read(gridStream);
00332
00333
00334 for (unsigned ntries=0; ntries<nFiles && g == 0; ++ntries)
00335 {
00336 gridStream.close();
00337 currentFileNum = (currentFileNum + 1U) % nFiles;
00338 gridStream.open(externalGridFiles[currentFileNum].c_str(),
00339 std::ios_base::in | std::ios_base::binary);
00340 if (!gridStream.is_open())
00341 throw cms::Exception("FFTJetBadConfig")
00342 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00343 " failed to open external grid file "
00344 << externalGridFiles[currentFileNum] << std::endl;
00345 g = fftjet::Grid2d<float>::read(gridStream);
00346 }
00347
00348 if (g)
00349 if (g->sum() > externalGridMaxEnergy)
00350 {
00351 delete g;
00352 g = 0;
00353 if (++nEnergyRejected >= maxFail)
00354 throw cms::Exception("FFTJetBadConfig")
00355 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00356 " too many grids in a row (" << nEnergyRejected
00357 << ") failed the maximum energy cut" << std::endl;
00358 }
00359 }
00360
00361 if (g)
00362 {
00363 add_Grid2d_data(energyFlow.get(), *g);
00364 delete g;
00365 }
00366 else
00367 {
00368
00369 throw cms::Exception("FFTJetBadConfig")
00370 << "ERROR in FFTJetPileupProcessor::mixExtraGrid():"
00371 " no valid grid records found" << std::endl;
00372 }
00373 }
00374
00375
00376
00377 void FFTJetPileupProcessor::beginJob()
00378 {
00379 }
00380
00381
00382
00383 void FFTJetPileupProcessor::endJob()
00384 {
00385 }
00386
00387
00388 void FFTJetPileupProcessor::loadFlatteningFactors(const edm::EventSetup& iSetup)
00389 {
00390 edm::ESHandle<FFTJetLookupTableSequence> h;
00391 StaticFFTJetLookupTableSequenceLoader::instance().load(
00392 iSetup, flatteningTableRecord, h);
00393 boost::shared_ptr<npstat::StorableMultivariateFunctor> f =
00394 (*h)[flatteningTableCategory][flatteningTableName];
00395
00396
00397 const unsigned nEta = energyFlow->nEta();
00398 etaFlatteningFactors.clear();
00399 etaFlatteningFactors.reserve(nEta);
00400 for (unsigned ieta=0; ieta<nEta; ++ieta)
00401 {
00402 const double eta = energyFlow->etaBinCenter(ieta);
00403 etaFlatteningFactors.push_back((*f)(&eta, 1U));
00404 }
00405 }
00406
00407
00408
00409 DEFINE_FWK_MODULE(FFTJetPileupProcessor);