CMS 3D CMS Logo

FFTJetProducer.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoJets/FFTJetProducers
4 // Class: FFTJetProducer
5 //
21 //
22 // Original Author: Igor Volobouev
23 // Created: Sun Jun 20 14:32:36 CDT 2010
24 //
25 //
26 
27 #ifndef RecoJets_FFTJetProducers_FFTJetProducer_h
28 #define RecoJets_FFTJetProducers_FFTJetProducer_h
29 
30 #include <memory>
31 
32 // FFTJet headers
33 #include "fftjet/AbsRecombinationAlg.hh"
34 #include "fftjet/AbsVectorRecombinationAlg.hh"
35 #include "fftjet/SparseClusteringTree.hh"
36 
37 // Additional FFTJet headers
38 #include "fftjet/VectorRecombinationAlgFactory.hh"
39 #include "fftjet/RecombinationAlgFactory.hh"
40 
49 
55 
57 
58 // framework include files
62 
64 
65 // local FFTJet-related definitions
69 
70 namespace fftjetcms {
72 }
73 
74 class CaloGeometry;
75 class CaloGeometryRecord;
76 class HcalTopology;
78 
79 //
80 // class declaration
81 //
83 public:
84  typedef fftjet::RecombinedJet<fftjetcms::VectorLike> RecoFFTJet;
85  typedef fftjet::SparseClusteringTree<fftjet::Peak, long> SparseTree;
86 
87  // Masks for the status bits. Do not add anything
88  // here -- higher bits (starting with 0x1000) will be
89  // used to indicate jet correction levels applied.
90  enum StatusBits {
91  RESOLUTION = 0xff,
96  };
97 
99 
100  explicit FFTJetProducer(const edm::ParameterSet&);
101  // Explicitly disable other ways to construct this object
102  FFTJetProducer() = delete;
103  FFTJetProducer(const FFTJetProducer&) = delete;
104  FFTJetProducer& operator=(const FFTJetProducer&) = delete;
105 
106  ~FFTJetProducer() override;
107 
108  // Parser for the resolution enum
110 
111 protected:
112  // Functions which should be overriden from the base
113  void beginJob() override;
114  void produce(edm::Event&, const edm::EventSetup&) override;
115  void endJob() override;
116 
117  // The following functions can be overriden by derived classes
118  // in order to adjust jet reconstruction algorithm behavior.
119 
120  // Override the following method in order to implement
121  // your own precluster selection strategy
122  virtual void selectPreclusters(const SparseTree& tree,
123  const fftjet::Functor1<bool, fftjet::Peak>& peakSelector,
124  std::vector<fftjet::Peak>* preclusters);
125 
126  // Precluster maker from GenJets (useful in calibration)
127  virtual void genJetPreclusters(const SparseTree& tree,
128  edm::Event&,
129  const edm::EventSetup&,
130  const fftjet::Functor1<bool, fftjet::Peak>& peakSelector,
131  std::vector<fftjet::Peak>* preclusters);
132 
133  // Override the following method (which by default does not do
134  // anything) in order to implement your own process-dependent
135  // assignment of membership functions to preclusters. This method
136  // will be called once per event, just before the main algorithm.
137  virtual void assignMembershipFunctions(std::vector<fftjet::Peak>* preclusters);
138 
139  // Parser for the peak selector
140  virtual std::unique_ptr<fftjet::Functor1<bool, fftjet::Peak> > parse_peakSelector(const edm::ParameterSet&);
141 
142  // Parser for the default jet membership function
143  virtual std::unique_ptr<fftjet::ScaleSpaceKernel> parse_jetMembershipFunction(const edm::ParameterSet&);
144 
145  // Parser for the background membership function
146  virtual std::unique_ptr<fftjetcms::AbsBgFunctor> parse_bgMembershipFunction(const edm::ParameterSet&);
147 
148  // Calculator for the recombination scale
149  virtual std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > parse_recoScaleCalcPeak(const edm::ParameterSet&);
150 
151  // Calculator for the recombination scale ratio
152  virtual std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > parse_recoScaleRatioCalcPeak(
153  const edm::ParameterSet&);
154 
155  // Calculator for the membership function factor
156  virtual std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > parse_memberFactorCalcPeak(const edm::ParameterSet&);
157 
158  // Similar calculators for the iterative algorithm
159  virtual std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > parse_recoScaleCalcJet(const edm::ParameterSet&);
160 
161  virtual std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > parse_recoScaleRatioCalcJet(const edm::ParameterSet&);
162 
163  virtual std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > parse_memberFactorCalcJet(const edm::ParameterSet&);
164 
165  // Calculator of the distance between jets which is used to make
166  // the decision about convergence of the iterative algorithm
167  virtual std::unique_ptr<fftjet::Functor2<double, RecoFFTJet, RecoFFTJet> > parse_jetDistanceCalc(
168  const edm::ParameterSet&);
169 
170  // Pile-up density calculator
171  virtual std::unique_ptr<fftjetcms::AbsPileupCalculator> parse_pileupDensityCalc(const edm::ParameterSet& ps);
172 
173  // The following function performs most of the precluster selection
174  // work in this module. You might want to reuse it if only a slight
175  // modification of the "selectPreclusters" method is desired.
176  void selectTreeNodes(const SparseTree& tree,
177  const fftjet::Functor1<bool, fftjet::Peak>& peakSelect,
178  std::vector<SparseTree::NodeId>* nodes);
179 
180 private:
181  typedef fftjet::AbsVectorRecombinationAlg<fftjetcms::VectorLike, fftjetcms::BgData> RecoAlg;
182  typedef fftjet::AbsRecombinationAlg<fftjetcms::Real, fftjetcms::VectorLike, fftjetcms::BgData> GridAlg;
183 
184  // Useful local utilities
185  template <class Real>
186  void loadSparseTreeData(const edm::Event&);
187 
188  void removeFakePreclusters();
189 
190  // The following methods do most of the work.
191  // The following function tells us if the grid was rebuilt.
192  bool loadEnergyFlow(const edm::Event& iEvent, std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow);
193  void buildGridAlg();
195  bool checkConvergence(const std::vector<RecoFFTJet>& previousIterResult, std::vector<RecoFFTJet>& thisIterResult);
198  void saveResults(edm::Event& iEvent, const edm::EventSetup&, unsigned nPreclustersFound);
199 
200  template <typename Jet>
202 
203  template <typename Jet>
204  void makeProduces(const std::string& alias, const std::string& tag);
205 
206  // The following function scans the pile-up density
207  // and fills the pile-up grid. Can be overriden if
208  // necessary.
210  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density);
211 
212  // Similar function for getting pile-up shape from the database
213  virtual void determinePileupDensityFromDB(const edm::Event& iEvent,
214  const edm::EventSetup& iSetup,
215  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density);
216 
217  // The following function builds the pile-up estimate
218  // for each jet
219  void determinePileup();
220 
221  // The following function returns the number of iterations
222  // performed. If this number equals to or less than "maxIterations"
223  // then the iterations have converged. If the number larger than
224  // "maxIterations" then the iterations failed to converge (note,
225  // however, that only "maxIterations" iterations would still be
226  // performed).
227  unsigned iterateJetReconstruction();
228 
229  // A function to set jet status bits
230  static void setJetStatusBit(RecoFFTJet* jet, int mask, bool value);
231 
232  //
233  // ----------member data ---------------------------
234  //
235 
236  // Local copy of the module configuration
238 
239  // Label for the tree produced by FFTJetPatRecoProducer
241 
242  // Are we going to use energy flow discretization grid as input
243  // to jet reconstruction?
245 
246  // Are we going to rebuild the energy flow discretization grid
247  // or to reuse the grid made by FFTJetPatRecoProducer?
248  const bool reuseExistingGrid;
249 
250  // Are we iterating?
251  const unsigned maxIterations;
252 
253  // Parameters which affect iteration convergence
254  const unsigned nJetsRequiredToConverge;
255  const double convergenceDistance;
256 
257  // Are we building assignments of collection members to jets?
258  const bool assignConstituents;
259 
260  // Are we resumming the constituents to determine jet 4-vectors?
261  // This might make sense if FFTJet is used in the crisp, gridded
262  // mode to determine jet areas, and vector recombination is desired.
263  const bool resumConstituents;
264 
265  // Are we going to subtract the pile-up? Note that
266  // pile-up subtraction does not modify eta and phi moments.
267  const bool calculatePileup;
268  const bool subtractPileup;
270 
271  // Label for the pile-up energy flow. Must be specified
272  // if the pile-up is subtracted.
274 
275  // Scale for the peak selection (if the scale is fixed)
276  const double fixedScale;
277 
278  // Minimum and maximum scale for searching stable configurations
279  const double minStableScale;
280  const double maxStableScale;
281 
282  // Stability "alpha"
283  const double stabilityAlpha;
284 
285  // Not sure at this point how to treat noise... For now, this is
286  // just a single configurable number...
287  const double noiseLevel;
288 
289  // Number of clusters requested (if the scale is adaptive)
290  const unsigned nClustersRequested;
291 
292  // Maximum eta for the grid-based algorithm
293  const double gridScanMaxEta;
294 
295  // Parameters related to the recombination algorithm
297  const bool isCrisp;
298  const double unlikelyBgWeight;
300 
301  // Label for the genJets used as seeds for jets
303 
304  // Maximum number of preclusters to use as jet seeds.
305  // This does not take into account the preclusters
306  // for which the value of the membership factor is 0.
307  const unsigned maxInitialPreclusters;
308 
309  // Resolution. The corresponding parameter value
310  // should be one of "fixed", "maximallyStable",
311  // "globallyAdaptive", "locallyAdaptive", or "fromGenJets".
313 
314  // Parameters related to the pileup shape stored
315  // in the database
320 
321  // Scales used
322  std::unique_ptr<std::vector<double> > iniScales;
323 
324  // The sparse clustering tree
326 
327  // Peak selector for the peaks already in the tree
328  std::unique_ptr<fftjet::Functor1<bool, fftjet::Peak> > peakSelector;
329 
330  // Recombination algorithms and related quantities
331  std::unique_ptr<RecoAlg> recoAlg;
332  std::unique_ptr<GridAlg> gridAlg;
333  std::unique_ptr<fftjet::ScaleSpaceKernel> jetMembershipFunction;
334  std::unique_ptr<fftjetcms::AbsBgFunctor> bgMembershipFunction;
335 
336  // Calculator for the recombination scale
337  std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > recoScaleCalcPeak;
338 
339  // Calculator for the recombination scale ratio
340  std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > recoScaleRatioCalcPeak;
341 
342  // Calculator for the membership function factor
343  std::unique_ptr<fftjet::Functor1<double, fftjet::Peak> > memberFactorCalcPeak;
344 
345  // Similar calculators for the iterative algorithm
346  std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > recoScaleCalcJet;
347  std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > recoScaleRatioCalcJet;
348  std::unique_ptr<fftjet::Functor1<double, RecoFFTJet> > memberFactorCalcJet;
349 
350  // Calculator for the jet distance used to estimate convergence
351  // of the iterative algorithm
352  std::unique_ptr<fftjet::Functor2<double, RecoFFTJet, RecoFFTJet> > jetDistanceCalc;
353 
354  // Vector of selected tree nodes
355  std::vector<SparseTree::NodeId> nodes;
356 
357  // Vector of selected preclusters
358  std::vector<fftjet::Peak> preclusters;
359 
360  // Vector of reconstructed jets (we will refill it in every event)
361  std::vector<RecoFFTJet> recoJets;
362 
363  // Cluster occupancy calculated as a function of level number
364  std::vector<unsigned> occupancy;
365 
366  // The thresholds obtained by the LOCALLY_ADAPTIVE method
367  std::vector<double> thresholds;
368 
369  // Minimum, maximum and used level calculated by some algorithms
371 
372  // Unclustered/unused energy produced during recombination
374  double unused;
375 
376  // Quantities defined below are used in the iterative mode only
377  std::vector<fftjet::Peak> iterPreclusters;
378  std::vector<RecoFFTJet> iterJets;
380 
381  // Vectors of constituents
382  std::vector<std::vector<reco::CandidatePtr> > constituents;
383 
384  // Vector of pile-up. We will subtract it from the
385  // 4-vectors of reconstructed jets.
386  std::vector<fftjetcms::VectorLike> pileup;
387 
388  // The pile-up transverse energy density discretization grid.
389  // Note that this is _density_, not energy. To get energy,
390  // multiply by cell area.
391  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> > pileupEnergyFlow;
392 
393  // The functor that calculates the pile-up density
394  std::unique_ptr<fftjetcms::AbsPileupCalculator> pileupDensityCalc;
395 
396  // Memory buffers related to pile-up subtraction
397  std::vector<fftjet::AbsKernel2d*> memFcns2dVec;
398  std::vector<double> doubleBuf;
399  std::vector<unsigned> cellCountsVec;
400 
401  // Tokens for data access
406 
409 };
410 
411 #endif // RecoJets_FFTJetProducers_FFTJetProducer_h
std::vector< std::vector< reco::CandidatePtr > > constituents
void produce(edm::Event &, const edm::EventSetup &) override
virtual std::unique_ptr< fftjet::ScaleSpaceKernel > parse_jetMembershipFunction(const edm::ParameterSet &)
std::unique_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > jetDistanceCalc
unsigned iterationsPerformed
const bool useGriddedAlgorithm
std::vector< fftjetcms::VectorLike > pileup
static Resolution parse_resolution(const std::string &name)
const bool resumConstituents
void selectTreeNodes(const SparseTree &tree, const fftjet::Functor1< bool, fftjet::Peak > &peakSelect, std::vector< SparseTree::NodeId > *nodes)
const double gridScanMaxEta
void saveResults(edm::Event &iEvent, const edm::EventSetup &, unsigned nPreclustersFound)
virtual void assignMembershipFunctions(std::vector< fftjet::Peak > *preclusters)
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > memberFactorCalcPeak
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleCalcJet
const std::string recombinationAlgorithm
const double fixedScale
const bool subtractPileupAs4Vec
virtual std::unique_ptr< fftjetcms::AbsBgFunctor > parse_bgMembershipFunction(const edm::ParameterSet &)
const double unlikelyBgWeight
void makeProduces(const std::string &alias, const std::string &tag)
const bool subtractPileup
const unsigned nJetsRequiredToConverge
Resolution resolution
std::vector< double > doubleBuf
std::vector< fftjet::AbsKernel2d * > memFcns2dVec
std::vector< RecoFFTJet > recoJets
std::unique_ptr< GridAlg > gridAlg
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geometry_token_
const bool isCrisp
const double recombinationDataCutoff
fftjet::AbsVectorRecombinationAlg< fftjetcms::VectorLike, fftjetcms::BgData > RecoAlg
std::vector< unsigned > occupancy
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
const double minStableScale
const double maxStableScale
std::vector< fftjet::Peak > preclusters
constexpr uint32_t mask
Definition: gpuClustering.h:24
FFTJetProducer & operator=(const FFTJetProducer &)=delete
void determineVectorConstituents()
FFTJetProducer()=delete
edm::EDGetTokenT< reco::PattRecoTree< fftjetcms::Real, reco::PattRecoPeak< fftjetcms::Real > > > input_recotree_token_
void loadSparseTreeData(const edm::Event &)
const bool reuseExistingGrid
const unsigned nClustersRequested
virtual std::unique_ptr< fftjetcms::AbsPileupCalculator > parse_pileupDensityCalc(const edm::ParameterSet &ps)
std::unique_ptr< std::vector< double > > iniScales
int iEvent
Definition: GenABIO.cc:224
fftjet::AbsRecombinationAlg< fftjetcms::Real, fftjetcms::VectorLike, fftjetcms::BgData > GridAlg
std::unique_ptr< fftjetcms::AbsBgFunctor > bgMembershipFunction
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleCalcJet(const edm::ParameterSet &)
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > memberFactorCalcJet
math::XYZTLorentzVector VectorLike
unsigned usedLevel
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleRatioCalcPeak(const edm::ParameterSet &)
void endJob() override
const double noiseLevel
static void setJetStatusBit(RecoFFTJet *jet, int mask, bool value)
std::vector< RecoFFTJet > iterJets
virtual void determinePileupDensityFromDB(const edm::Event &iEvent, const edm::EventSetup &iSetup, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
std::vector< unsigned > cellCountsVec
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleRatioCalcJet
std::unique_ptr< fftjet::ScaleSpaceKernel > jetMembershipFunction
Definition: value.py:1
std::string pileupTableRecord
std::string pileupTableName
std::vector< double > thresholds
unsigned iterateJetReconstruction()
const bool calculatePileup
const bool assignConstituents
bool checkConvergence(const std::vector< RecoFFTJet > &previousIterResult, std::vector< RecoFFTJet > &thisIterResult)
A grid filled with discretized energy flow.
const edm::InputTag treeLabel
virtual void selectPreclusters(const SparseTree &tree, const fftjet::Functor1< bool, fftjet::Peak > &peakSelector, std::vector< fftjet::Peak > *preclusters)
const edm::InputTag genJetsLabel
const edm::ParameterSet myConfiguration
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > input_energyflow_token_
fftjet::RecombinedJet< fftjetcms::VectorLike > RecoFFTJet
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_recoScaleRatioCalcJet(const edm::ParameterSet &)
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > topology_token_
edm::EDGetTokenT< reco::FFTJetPileupSummary > input_pusummary_token_
const double stabilityAlpha
const unsigned maxInitialPreclusters
std::string pileupTableCategory
SparseTree sparseTree
std::unique_ptr< fftjetcms::AbsPileupCalculator > pileupDensityCalc
std::vector< SparseTree::NodeId > nodes
void determineGriddedConstituents()
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
const double convergenceDistance
edm::EDGetTokenT< std::vector< reco::FFTAnyJet< reco::GenJet > > > input_genjet_token_
void writeJets(edm::Event &iEvent, const edm::EventSetup &)
~FFTJetProducer() override
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleRatioCalcPeak
const edm::InputTag pileupLabel
void beginJob() override
virtual std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > parse_peakSelector(const edm::ParameterSet &)
Definition: tree.py:1
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > pileupEnergyFlow
virtual void genJetPreclusters(const SparseTree &tree, edm::Event &, const edm::EventSetup &, const fftjet::Functor1< bool, fftjet::Peak > &peakSelector, std::vector< fftjet::Peak > *preclusters)
bool loadEnergyFlow(const edm::Event &iEvent, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &flow)
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleCalcPeak
std::unique_ptr< RecoAlg > recoAlg
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_memberFactorCalcPeak(const edm::ParameterSet &)
fftjetcms::VectorLike unclustered
virtual std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > parse_recoScaleCalcPeak(const edm::ParameterSet &)
std::vector< fftjet::Peak > iterPreclusters
virtual void determinePileupDensityFromConfig(const edm::Event &iEvent, std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > &density)
void removeFakePreclusters()
const unsigned maxIterations
void prepareRecombinationScales()
virtual std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > parse_memberFactorCalcJet(const edm::ParameterSet &)
virtual std::unique_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > parse_jetDistanceCalc(const edm::ParameterSet &)