CMS 3D CMS Logo

FFTJetPileupEstimator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoJets/FFTJetProducers
4 // Class: FFTJetPileupEstimator
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Wed Apr 20 13:52:23 CDT 2011
16 //
17 //
18 
19 #include <cmath>
20 
21 // Framework include files
26 
27 // Data formats
30 // #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
33 
35 
36 // Loader for the lookup tables
38 
39 #define init_param(type, varname) varname(ps.getParameter<type>(#varname))
40 
41 using namespace fftjetcms;
42 
43 //
44 // class declaration
45 //
47 public:
49  FFTJetPileupEstimator() = delete;
51  FFTJetPileupEstimator& operator=(const FFTJetPileupEstimator&) = delete;
52  ~FFTJetPileupEstimator() override;
53 
54 protected:
55  // methods
56  void beginJob() override;
57  void produce(edm::Event&, const edm::EventSetup&) override;
58  void endJob() override;
59 
60 private:
61  std::unique_ptr<reco::FFTJetPileupSummary> calibrateFromConfig(double uncalibrated) const;
62 
63  std::unique_ptr<reco::FFTJetPileupSummary> calibrateFromDB(double uncalibrated, const edm::EventSetup& iSetup) const;
64 
65  template <class Ptr>
66  inline void checkConfig(const Ptr& ptr, const char* message) {
67  if (ptr.get() == nullptr)
68  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
69  }
70 
73 
75  double cdfvalue;
77  unsigned filterNumber;
78  std::vector<double> uncertaintyZones;
79  std::unique_ptr<fftjet::Functor1<double, double> > calibrationCurve;
80  std::unique_ptr<fftjet::Functor1<double, double> > uncertaintyCurve;
81 
82  // Alternative method to calibrate the pileup.
83  // We will fetch three lookup tables from the database:
84  //
85  // 1. calibration curve table
86  //
87  // 2. uncertainty curve table
88  //
89  // 3. the table that will lookup the uncertainty code
90  // given the uncalibrated pileup value
91  //
92  // It will be assumed that all these tables will have
93  // the same record and category but different names.
94  //
101 };
102 
103 //
104 // constructors and destructor
105 //
109  init_param(double, cdfvalue),
110  init_param(double, ptToDensityFactor),
111  init_param(unsigned, filterNumber),
120  checkConfig(calibrationCurve, "bad calibration curve definition");
121 
123  checkConfig(uncertaintyCurve, "bad uncertainty curve definition");
124 
125  inputToken = consumes<reco::DiscretizedEnergyFlow>(inputLabel);
126 
127  produces<reco::FFTJetPileupSummary>(outputLabel);
128 }
129 
131 
132 //
133 // member functions
134 //
135 
136 // ------------ method called to for each event ------------
139  iEvent.getByToken(inputToken, input);
140 
142  const unsigned nScales = h.nEtaBins();
143  const unsigned nCdfvalues = h.nPhiBins();
144 
145  const unsigned fixedCdfvalueBin = static_cast<unsigned>(std::floor(cdfvalue * nCdfvalues));
146  if (fixedCdfvalueBin >= nCdfvalues) {
147  throw cms::Exception("FFTJetBadConfig") << "Bad cdf value" << std::endl;
148  }
149  if (filterNumber >= nScales) {
150  throw cms::Exception("FFTJetBadConfig") << "Bad filter number" << std::endl;
151  }
152 
153  // Simple fixed-point pile-up estimate
154  const double curve = h.data()[filterNumber * nCdfvalues + fixedCdfvalueBin];
155 
156  std::unique_ptr<reco::FFTJetPileupSummary> summary;
157  if (loadCalibFromDB)
158  summary = calibrateFromDB(curve, iSetup);
159  else
160  summary = calibrateFromConfig(curve);
162 }
163 
165 
167 
168 std::unique_ptr<reco::FFTJetPileupSummary> FFTJetPileupEstimator::calibrateFromConfig(const double curve) const {
169  const double pileupRho = ptToDensityFactor * (*calibrationCurve)(curve);
170  const double rhoUncert = ptToDensityFactor * (*uncertaintyCurve)(curve);
171 
172  // Determine the uncertainty zone of the estimate. The "curve"
173  // has to be above or equal to uncertaintyZones[i] but below
174  // uncertaintyZones[i + 1] (the second condition is also satisfied
175  // by i == uncertaintyZones.size() - 1). Of course, it is assumed
176  // that the vector of zones is configured appropriately -- the zone
177  // boundaries must be presented in the increasing order.
178  int uncertaintyCode = -1;
179  if (!uncertaintyZones.empty()) {
180  const unsigned nZones = uncertaintyZones.size();
181  for (unsigned i = 0; i < nZones; ++i)
182  if (curve >= uncertaintyZones[i]) {
183  if (i == nZones - 1U) {
184  uncertaintyCode = i;
185  break;
186  } else if (curve < uncertaintyZones[i + 1]) {
187  uncertaintyCode = i;
188  break;
189  }
190  }
191  }
192 
193  return std::make_unique<reco::FFTJetPileupSummary>(curve, pileupRho, rhoUncert, uncertaintyCode);
194 }
195 
196 std::unique_ptr<reco::FFTJetPileupSummary> FFTJetPileupEstimator::calibrateFromDB(const double curve,
197  const edm::EventSetup& iSetup) const {
200  std::shared_ptr<npstat::StorableMultivariateFunctor> uz = (*h)[calibTableCategory][uncertaintyZonesName];
201  std::shared_ptr<npstat::StorableMultivariateFunctor> cc = (*h)[calibTableCategory][calibrationCurveName];
202  std::shared_ptr<npstat::StorableMultivariateFunctor> uc = (*h)[calibTableCategory][uncertaintyCurveName];
203 
204  const double pileupRho = ptToDensityFactor * (*cc)(&curve, 1U);
205  const double rhoUncert = ptToDensityFactor * (*uc)(&curve, 1U);
206  const int uncertaintyCode = round((*uz)(&curve, 1U));
207 
208  return std::make_unique<reco::FFTJetPileupSummary>(curve, pileupRho, rhoUncert, uncertaintyCode);
209 }
210 
211 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::Event &, const edm::EventSetup &) override
FFTJetPileupEstimator()=delete
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::unique_ptr< reco::FFTJetPileupSummary > calibrateFromConfig(double uncalibrated) const
void beginJob()
Definition: Breakpoints.cc:14
static std::string const input
Definition: EdmProvDump.cc:47
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > inputToken
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< fftjet::Functor1< double, double > > calibrationCurve
std::vector< double > uncertaintyZones
std::unique_ptr< fftjet::Functor1< double, double > > fftjet_Function_parser(const edm::ParameterSet &ps)
std::unique_ptr< reco::FFTJetPileupSummary > calibrateFromDB(double uncalibrated, const edm::EventSetup &iSetup) const
HLT enums.
std::unique_ptr< fftjet::Functor1< double, double > > uncertaintyCurve
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
static const Mapper & instance()
#define init_param(type, varname)
def move(src, dest)
Definition: eostools.py:511
void checkConfig(const Ptr &ptr, const char *message)