CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
48 public:
51 
52 protected:
53  // methods
54  void beginJob() override;
55  void produce(edm::Event&, const edm::EventSetup&) override;
56  void endJob() override;
57 
58 private:
62 
63  std::auto_ptr<reco::FFTJetPileupSummary> calibrateFromConfig(
64  double uncalibrated) const;
65 
66  std::auto_ptr<reco::FFTJetPileupSummary> calibrateFromDB(
67  double uncalibrated, const edm::EventSetup& iSetup) const;
68 
69  template<class Ptr>
70  inline void checkConfig(const Ptr& ptr, const char* message)
71  {
72  if (ptr.get() == NULL)
73  throw cms::Exception("FFTJetBadConfig") << message << std::endl;
74  }
75 
78  double cdfvalue;
80  unsigned filterNumber;
81  std::vector<double> uncertaintyZones;
82  std::auto_ptr<fftjet::Functor1<double,double> > calibrationCurve;
83  std::auto_ptr<fftjet::Functor1<double,double> > uncertaintyCurve;
84 
85  // Alternative method to calibrate the pileup.
86  // We will fetch three lookup tables from the database:
87  //
88  // 1. calibration curve table
89  //
90  // 2. uncertainty curve table
91  //
92  // 3. the table that will lookup the uncertainty code
93  // given the uncalibrated pileup value
94  //
95  // It will be assumed that all these tables will have
96  // the same record and category but different names.
97  //
104 };
105 
106 //
107 // constructors and destructor
108 //
110  : init_param(edm::InputTag, inputLabel),
111  init_param(std::string, outputLabel),
112  init_param(double, cdfvalue),
113  init_param(double, ptToDensityFactor),
114  init_param(unsigned, filterNumber),
115  init_param(std::vector<double>, uncertaintyZones),
116  init_param(std::string, calibTableRecord),
117  init_param(std::string, calibTableCategory),
118  init_param(std::string, uncertaintyZonesName),
119  init_param(std::string, calibrationCurveName),
120  init_param(std::string, uncertaintyCurveName),
121  init_param(bool, loadCalibFromDB)
122 {
124  ps.getParameter<edm::ParameterSet>("calibrationCurve"));
125  checkConfig(calibrationCurve, "bad calibration curve definition");
126 
128  ps.getParameter<edm::ParameterSet>("uncertaintyCurve"));
129  checkConfig(uncertaintyCurve, "bad uncertainty curve definition");
130 
131  produces<reco::FFTJetPileupSummary>(outputLabel);
132 }
133 
134 
136 {
137 }
138 
139 //
140 // member functions
141 //
142 
143 // ------------ method called to for each event ------------
145  const edm::EventSetup& iSetup)
146 {
148  iEvent.getByLabel(inputLabel, input);
149 
150  const reco::DiscretizedEnergyFlow& h(*input);
151  const unsigned nScales = h.nEtaBins();
152  const unsigned nCdfvalues = h.nPhiBins();
153 
154  const unsigned fixedCdfvalueBin = static_cast<unsigned>(
155  std::floor(cdfvalue*nCdfvalues));
156  if (fixedCdfvalueBin >= nCdfvalues)
157  {
158  throw cms::Exception("FFTJetBadConfig")
159  << "Bad cdf value" << std::endl;
160  }
161  if (filterNumber >= nScales)
162  {
163  throw cms::Exception("FFTJetBadConfig")
164  << "Bad filter number" << std::endl;
165  }
166 
167  // Simple fixed-point pile-up estimate
168  const double curve = h.data()[filterNumber*nCdfvalues + fixedCdfvalueBin];
169 
170  std::auto_ptr<reco::FFTJetPileupSummary> summary;
171  if (loadCalibFromDB)
172  summary = calibrateFromDB(curve, iSetup);
173  else
174  summary = calibrateFromConfig(curve);
175  iEvent.put(summary, outputLabel);
176 }
177 
178 
180 {
181 }
182 
183 
185 {
186 }
187 
188 
189 std::auto_ptr<reco::FFTJetPileupSummary>
191 {
192  const double pileupRho = ptToDensityFactor*(*calibrationCurve)(curve);
193  const double rhoUncert = ptToDensityFactor*(*uncertaintyCurve)(curve);
194 
195  // Determine the uncertainty zone of the estimate. The "curve"
196  // has to be above or equal to uncertaintyZones[i] but below
197  // uncertaintyZones[i + 1] (the second condition is also satisfied
198  // by i == uncertaintyZones.size() - 1). Of course, it is assumed
199  // that the vector of zones is configured appropriately -- the zone
200  // boundaries must be presented in the increasing order.
201  int uncertaintyCode = -1;
202  if (!uncertaintyZones.empty())
203  {
204  const unsigned nZones = uncertaintyZones.size();
205  for (unsigned i = 0; i < nZones; ++i)
206  if (curve >= uncertaintyZones[i])
207  {
208  if (i == nZones - 1U)
209  {
210  uncertaintyCode = i;
211  break;
212  }
213  else if (curve < uncertaintyZones[i + 1])
214  {
215  uncertaintyCode = i;
216  break;
217  }
218  }
219  }
220 
221  return std::auto_ptr<reco::FFTJetPileupSummary>(
222  new reco::FFTJetPileupSummary(curve, pileupRho,
223  rhoUncert, uncertaintyCode));
224 }
225 
226 
227 std::auto_ptr<reco::FFTJetPileupSummary>
229  const double curve, const edm::EventSetup& iSetup) const
230 {
233  iSetup, calibTableRecord, h);
234  boost::shared_ptr<npstat::StorableMultivariateFunctor> uz =
236  boost::shared_ptr<npstat::StorableMultivariateFunctor> cc =
238  boost::shared_ptr<npstat::StorableMultivariateFunctor> uc =
240 
241  const double pileupRho = ptToDensityFactor*(*cc)(&curve, 1U);
242  const double rhoUncert = ptToDensityFactor*(*uc)(&curve, 1U);
243  const int uncertaintyCode = round((*uz)(&curve, 1U));
244 
245  return std::auto_ptr<reco::FFTJetPileupSummary>(
246  new reco::FFTJetPileupSummary(curve, pileupRho,
247  rhoUncert, uncertaintyCode));
248 }
249 
250 
251 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
#define NULL
Definition: scimark2.h:8
Summary info for pile-up determined by Gaussian filtering.
std::auto_ptr< fftjet::Functor1< double, double > > calibrationCurve
void beginJob()
Definition: Breakpoints.cc:15
static std::string const input
Definition: EdmProvDump.cc:44
const double * data() const
int iEvent
Definition: GenABIO.cc:243
std::auto_ptr< fftjet::Functor1< double, double > > fftjet_Function_parser(const edm::ParameterSet &ps)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
std::vector< double > uncertaintyZones
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
std::auto_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
std::auto_ptr< reco::FFTJetPileupSummary > calibrateFromDB(double uncalibrated, const edm::EventSetup &iSetup) const
static const Mapper & instance()
#define init_param(type, varname)
void checkConfig(const Ptr &ptr, const char *message)
std::auto_ptr< reco::FFTJetPileupSummary > calibrateFromConfig(double uncalibrated) const