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 
80  double cdfvalue;
82  unsigned filterNumber;
83  std::vector<double> uncertaintyZones;
84  std::auto_ptr<fftjet::Functor1<double,double> > calibrationCurve;
85  std::auto_ptr<fftjet::Functor1<double,double> > uncertaintyCurve;
86 
87  // Alternative method to calibrate the pileup.
88  // We will fetch three lookup tables from the database:
89  //
90  // 1. calibration curve table
91  //
92  // 2. uncertainty curve table
93  //
94  // 3. the table that will lookup the uncertainty code
95  // given the uncalibrated pileup value
96  //
97  // It will be assumed that all these tables will have
98  // the same record and category but different names.
99  //
106 };
107 
108 //
109 // constructors and destructor
110 //
112  : init_param(edm::InputTag, inputLabel),
113  init_param(std::string, outputLabel),
114  init_param(double, cdfvalue),
115  init_param(double, ptToDensityFactor),
116  init_param(unsigned, filterNumber),
117  init_param(std::vector<double>, uncertaintyZones),
118  init_param(std::string, calibTableRecord),
119  init_param(std::string, calibTableCategory),
120  init_param(std::string, uncertaintyZonesName),
121  init_param(std::string, calibrationCurveName),
122  init_param(std::string, uncertaintyCurveName),
123  init_param(bool, loadCalibFromDB)
124 {
126  ps.getParameter<edm::ParameterSet>("calibrationCurve"));
127  checkConfig(calibrationCurve, "bad calibration curve definition");
128 
130  ps.getParameter<edm::ParameterSet>("uncertaintyCurve"));
131  checkConfig(uncertaintyCurve, "bad uncertainty curve definition");
132 
133  inputToken = consumes<reco::DiscretizedEnergyFlow>(inputLabel);
134 
135  produces<reco::FFTJetPileupSummary>(outputLabel);
136 }
137 
138 
140 {
141 }
142 
143 //
144 // member functions
145 //
146 
147 // ------------ method called to for each event ------------
149  const edm::EventSetup& iSetup)
150 {
152  iEvent.getByToken(inputToken, input);
153 
154  const reco::DiscretizedEnergyFlow& h(*input);
155  const unsigned nScales = h.nEtaBins();
156  const unsigned nCdfvalues = h.nPhiBins();
157 
158  const unsigned fixedCdfvalueBin = static_cast<unsigned>(
159  std::floor(cdfvalue*nCdfvalues));
160  if (fixedCdfvalueBin >= nCdfvalues)
161  {
162  throw cms::Exception("FFTJetBadConfig")
163  << "Bad cdf value" << std::endl;
164  }
165  if (filterNumber >= nScales)
166  {
167  throw cms::Exception("FFTJetBadConfig")
168  << "Bad filter number" << std::endl;
169  }
170 
171  // Simple fixed-point pile-up estimate
172  const double curve = h.data()[filterNumber*nCdfvalues + fixedCdfvalueBin];
173 
174  std::auto_ptr<reco::FFTJetPileupSummary> summary;
175  if (loadCalibFromDB)
176  summary = calibrateFromDB(curve, iSetup);
177  else
178  summary = calibrateFromConfig(curve);
179  iEvent.put(summary, outputLabel);
180 }
181 
182 
184 {
185 }
186 
187 
189 {
190 }
191 
192 
193 std::auto_ptr<reco::FFTJetPileupSummary>
195 {
196  const double pileupRho = ptToDensityFactor*(*calibrationCurve)(curve);
197  const double rhoUncert = ptToDensityFactor*(*uncertaintyCurve)(curve);
198 
199  // Determine the uncertainty zone of the estimate. The "curve"
200  // has to be above or equal to uncertaintyZones[i] but below
201  // uncertaintyZones[i + 1] (the second condition is also satisfied
202  // by i == uncertaintyZones.size() - 1). Of course, it is assumed
203  // that the vector of zones is configured appropriately -- the zone
204  // boundaries must be presented in the increasing order.
205  int uncertaintyCode = -1;
206  if (!uncertaintyZones.empty())
207  {
208  const unsigned nZones = uncertaintyZones.size();
209  for (unsigned i = 0; i < nZones; ++i)
210  if (curve >= uncertaintyZones[i])
211  {
212  if (i == nZones - 1U)
213  {
214  uncertaintyCode = i;
215  break;
216  }
217  else if (curve < uncertaintyZones[i + 1])
218  {
219  uncertaintyCode = i;
220  break;
221  }
222  }
223  }
224 
225  return std::auto_ptr<reco::FFTJetPileupSummary>(
226  new reco::FFTJetPileupSummary(curve, pileupRho,
227  rhoUncert, uncertaintyCode));
228 }
229 
230 
231 std::auto_ptr<reco::FFTJetPileupSummary>
233  const double curve, const edm::EventSetup& iSetup) const
234 {
237  iSetup, calibTableRecord, h);
238  boost::shared_ptr<npstat::StorableMultivariateFunctor> uz =
240  boost::shared_ptr<npstat::StorableMultivariateFunctor> cc =
242  boost::shared_ptr<npstat::StorableMultivariateFunctor> uc =
244 
245  const double pileupRho = ptToDensityFactor*(*cc)(&curve, 1U);
246  const double rhoUncert = ptToDensityFactor*(*uc)(&curve, 1U);
247  const int uncertaintyCode = round((*uz)(&curve, 1U));
248 
249  return std::auto_ptr<reco::FFTJetPileupSummary>(
250  new reco::FFTJetPileupSummary(curve, pileupRho,
251  rhoUncert, uncertaintyCode));
252 }
253 
254 
255 //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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#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:43
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > inputToken
const double * data() const
int iEvent
Definition: GenABIO.cc:230
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:120
std::vector< double > uncertaintyZones
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< fftjet::Functor1< double, double > > uncertaintyCurve
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