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 
59 // framework include files
63 
65 
66 // local FFTJet-related definitions
70 
71 namespace fftjetcms {
73 }
74 
75 //
76 // class declaration
77 //
79 {
80 public:
81  typedef fftjet::RecombinedJet<fftjetcms::VectorLike> RecoFFTJet;
82  typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
83 
84  // Masks for the status bits. Do not add anything
85  // here -- higher bits (starting with 0x1000) will be
86  // used to indicate jet correction levels applied.
88  {
89  RESOLUTION = 0xff,
90  CONSTITUENTS_RESUMMED = 0x100,
91  PILEUP_CALCULATED = 0x200,
92  PILEUP_SUBTRACTED_4VEC = 0x400,
93  PILEUP_SUBTRACTED_PT = 0x800
94  };
95 
97  {
98  FIXED = 0,
102  FROM_GENJETS
103  };
104 
105  explicit FFTJetProducer(const edm::ParameterSet&);
106  ~FFTJetProducer() override;
107 
108  // Parser for the resolution enum
109  static Resolution parse_resolution(const std::string& name);
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(
123  const SparseTree& tree,
124  const fftjet::Functor1<bool,fftjet::Peak>& peakSelector,
125  std::vector<fftjet::Peak>* preclusters);
126 
127  // Precluster maker from GenJets (useful in calibration)
128  virtual void genJetPreclusters(
129  const SparseTree& tree,
130  edm::Event&, const edm::EventSetup&,
131  const fftjet::Functor1<bool,fftjet::Peak>& peakSelector,
132  std::vector<fftjet::Peak>* preclusters);
133 
134  // Override the following method (which by default does not do
135  // anything) in order to implement your own process-dependent
136  // assignment of membership functions to preclusters. This method
137  // will be called once per event, just before the main algorithm.
138  virtual void assignMembershipFunctions(
139  std::vector<fftjet::Peak>* preclusters);
140 
141  // Parser for the peak selector
142  virtual std::unique_ptr<fftjet::Functor1<bool,fftjet::Peak> >
143  parse_peakSelector(const edm::ParameterSet&);
144 
145  // Parser for the default jet membership function
146  virtual std::unique_ptr<fftjet::ScaleSpaceKernel>
147  parse_jetMembershipFunction(const edm::ParameterSet&);
148 
149  // Parser for the background membership function
150  virtual std::unique_ptr<fftjetcms::AbsBgFunctor>
151  parse_bgMembershipFunction(const edm::ParameterSet&);
152 
153  // Calculator for the recombination scale
154  virtual std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> >
155  parse_recoScaleCalcPeak(const edm::ParameterSet&);
156 
157  // Calculator for the recombination scale ratio
158  virtual std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> >
159  parse_recoScaleRatioCalcPeak(const edm::ParameterSet&);
160 
161  // Calculator for the membership function factor
162  virtual std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> >
163  parse_memberFactorCalcPeak(const edm::ParameterSet&);
164 
165  // Similar calculators for the iterative algorithm
166  virtual std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> >
167  parse_recoScaleCalcJet(const edm::ParameterSet&);
168 
169  virtual std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> >
170  parse_recoScaleRatioCalcJet(const edm::ParameterSet&);
171 
172  virtual std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> >
173  parse_memberFactorCalcJet(const edm::ParameterSet&);
174 
175  // Calculator of the distance between jets which is used to make
176  // the decision about convergence of the iterative algorithm
177  virtual std::unique_ptr<fftjet::Functor2<double,RecoFFTJet,RecoFFTJet> >
178  parse_jetDistanceCalc(const edm::ParameterSet&);
179 
180  // Pile-up density calculator
181  virtual std::unique_ptr<fftjetcms::AbsPileupCalculator>
182  parse_pileupDensityCalc(const edm::ParameterSet& ps);
183 
184  // The following function performs most of the precluster selection
185  // work in this module. You might want to reuse it if only a slight
186  // modification of the "selectPreclusters" method is desired.
187  void selectTreeNodes(const SparseTree& tree,
188  const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
189  std::vector<SparseTree::NodeId>* nodes);
190 private:
191  typedef fftjet::AbsVectorRecombinationAlg<
193  typedef fftjet::AbsRecombinationAlg<
195 
196  // Explicitly disable other ways to construct this object
197  FFTJetProducer() = delete;
198  FFTJetProducer(const FFTJetProducer&) = delete;
199  FFTJetProducer& operator=(const FFTJetProducer&) = delete;
200 
201  // Useful local utilities
202  template<class Real>
203  void loadSparseTreeData(const edm::Event&);
204 
205  void removeFakePreclusters();
206 
207  // The following methods do most of the work.
208  // The following function tells us if the grid was rebuilt.
209  bool loadEnergyFlow(
210  const edm::Event& iEvent,
211  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow);
212  void buildGridAlg();
213  void prepareRecombinationScales();
214  bool checkConvergence(const std::vector<RecoFFTJet>& previousIterResult,
215  std::vector<RecoFFTJet>& thisIterResult);
216  void determineGriddedConstituents();
217  void determineVectorConstituents();
218  void saveResults(edm::Event& iEvent, const edm::EventSetup&,
219  unsigned nPreclustersFound);
220 
221  template <typename Jet>
222  void writeJets(edm::Event& iEvent, const edm::EventSetup&);
223 
224  template <typename Jet>
225  void makeProduces(const std::string& alias, const std::string& tag);
226 
227  // The following function scans the pile-up density
228  // and fills the pile-up grid. Can be overriden if
229  // necessary.
230  virtual void determinePileupDensityFromConfig(
231  const edm::Event& iEvent,
232  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density);
233 
234  // Similar function for getting pile-up shape from the database
235  virtual void determinePileupDensityFromDB(
236  const edm::Event& iEvent, const edm::EventSetup& iSetup,
237  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> >& density);
238 
239  // The following function builds the pile-up estimate
240  // for each jet
241  void determinePileup();
242 
243  // The following function returns the number of iterations
244  // performed. If this number equals to or less than "maxIterations"
245  // then the iterations have converged. If the number larger than
246  // "maxIterations" then the iterations failed to converge (note,
247  // however, that only "maxIterations" iterations would still be
248  // performed).
249  unsigned iterateJetReconstruction();
250 
251  // A function to set jet status bits
252  static void setJetStatusBit(RecoFFTJet* jet, int mask, bool value);
253 
254  //
255  // ----------member data ---------------------------
256  //
257 
258  // Local copy of the module configuration
260 
261  // Label for the tree produced by FFTJetPatRecoProducer
263 
264  // Are we going to use energy flow discretization grid as input
265  // to jet reconstruction?
267 
268  // Are we going to rebuild the energy flow discretization grid
269  // or to reuse the grid made by FFTJetPatRecoProducer?
270  const bool reuseExistingGrid;
271 
272  // Are we iterating?
273  const unsigned maxIterations;
274 
275  // Parameters which affect iteration convergence
276  const unsigned nJetsRequiredToConverge;
277  const double convergenceDistance;
278 
279  // Are we building assignments of collection members to jets?
280  const bool assignConstituents;
281 
282  // Are we resumming the constituents to determine jet 4-vectors?
283  // This might make sense if FFTJet is used in the crisp, gridded
284  // mode to determine jet areas, and vector recombination is desired.
285  const bool resumConstituents;
286 
287  // Are we going to subtract the pile-up? Note that
288  // pile-up subtraction does not modify eta and phi moments.
289  const bool calculatePileup;
290  const bool subtractPileup;
292 
293  // Label for the pile-up energy flow. Must be specified
294  // if the pile-up is subtracted.
296 
297  // Scale for the peak selection (if the scale is fixed)
298  const double fixedScale;
299 
300  // Minimum and maximum scale for searching stable configurations
301  const double minStableScale;
302  const double maxStableScale;
303 
304  // Stability "alpha"
305  const double stabilityAlpha;
306 
307  // Not sure at this point how to treat noise... For now, this is
308  // just a single configurable number...
309  const double noiseLevel;
310 
311  // Number of clusters requested (if the scale is adaptive)
312  const unsigned nClustersRequested;
313 
314  // Maximum eta for the grid-based algorithm
315  const double gridScanMaxEta;
316 
317  // Parameters related to the recombination algorithm
319  const bool isCrisp;
320  const double unlikelyBgWeight;
322 
323  // Label for the genJets used as seeds for jets
325 
326  // Maximum number of preclusters to use as jet seeds.
327  // This does not take into account the preclusters
328  // for which the value of the membership factor is 0.
329  const unsigned maxInitialPreclusters;
330 
331  // Resolution. The corresponding parameter value
332  // should be one of "fixed", "maximallyStable",
333  // "globallyAdaptive", "locallyAdaptive", or "fromGenJets".
335 
336  // Parameters related to the pileup shape stored
337  // in the database
342 
343  // Scales used
344  std::unique_ptr<std::vector<double> > iniScales;
345 
346  // The sparse clustering tree
347  SparseTree sparseTree;
348 
349  // Peak selector for the peaks already in the tree
350  std::unique_ptr<fftjet::Functor1<bool,fftjet::Peak> > peakSelector;
351 
352  // Recombination algorithms and related quantities
353  std::unique_ptr<RecoAlg> recoAlg;
354  std::unique_ptr<GridAlg> gridAlg;
355  std::unique_ptr<fftjet::ScaleSpaceKernel> jetMembershipFunction;
356  std::unique_ptr<fftjetcms::AbsBgFunctor> bgMembershipFunction;
357 
358  // Calculator for the recombination scale
359  std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> > recoScaleCalcPeak;
360 
361  // Calculator for the recombination scale ratio
362  std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> >
364 
365  // Calculator for the membership function factor
366  std::unique_ptr<fftjet::Functor1<double,fftjet::Peak> > memberFactorCalcPeak;
367 
368  // Similar calculators for the iterative algorithm
369  std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> > recoScaleCalcJet;
370  std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> > recoScaleRatioCalcJet;
371  std::unique_ptr<fftjet::Functor1<double,RecoFFTJet> > memberFactorCalcJet;
372 
373  // Calculator for the jet distance used to estimate convergence
374  // of the iterative algorithm
375  std::unique_ptr<fftjet::Functor2<double,RecoFFTJet,RecoFFTJet> >
377 
378  // Vector of selected tree nodes
379  std::vector<SparseTree::NodeId> nodes;
380 
381  // Vector of selected preclusters
382  std::vector<fftjet::Peak> preclusters;
383 
384  // Vector of reconstructed jets (we will refill it in every event)
385  std::vector<RecoFFTJet> recoJets;
386 
387  // Cluster occupancy calculated as a function of level number
388  std::vector<unsigned> occupancy;
389 
390  // The thresholds obtained by the LOCALLY_ADAPTIVE method
391  std::vector<double> thresholds;
392 
393  // Minimum, maximum and used level calculated by some algorithms
394  unsigned minLevel, maxLevel, usedLevel;
395 
396  // Unclustered/unused energy produced during recombination
397  fftjetcms::VectorLike unclustered;
398  double unused;
399 
400  // Quantities defined below are used in the iterative mode only
401  std::vector<fftjet::Peak> iterPreclusters;
402  std::vector<RecoFFTJet> iterJets;
404 
405  // Vectors of constituents
406  std::vector<std::vector<reco::CandidatePtr> > constituents;
407 
408  // Vector of pile-up. We will subtract it from the
409  // 4-vectors of reconstructed jets.
410  std::vector<fftjetcms::VectorLike> pileup;
411 
412  // The pile-up transverse energy density discretization grid.
413  // Note that this is _density_, not energy. To get energy,
414  // multiply by cell area.
415  std::unique_ptr<fftjet::Grid2d<fftjetcms::Real> > pileupEnergyFlow;
416 
417  // The functor that calculates the pile-up density
418  std::unique_ptr<fftjetcms::AbsPileupCalculator> pileupDensityCalc;
419 
420  // Memory buffers related to pile-up subtraction
421  std::vector<fftjet::AbsKernel2d*> memFcns2dVec;
422  std::vector<double> doubleBuf;
423  std::vector<unsigned> cellCountsVec;
424 
425  // Tokens for data access
430 };
431 
432 #endif // RecoJets_FFTJetProducers_FFTJetProducer_h
std::vector< std::vector< reco::CandidatePtr > > constituents
unsigned iterationsPerformed
const bool useGriddedAlgorithm
std::vector< fftjetcms::VectorLike > pileup
const bool resumConstituents
const double gridScanMaxEta
double BgData
const std::string recombinationAlgorithm
const double fixedScale
const bool subtractPileupAs4Vec
const double unlikelyBgWeight
std::unique_ptr< fftjet::Functor1< bool, fftjet::Peak > > peakSelector
const bool subtractPileup
const unsigned nJetsRequiredToConverge
Resolution resolution
edm::EDGetTokenT< reco::PattRecoTree< fftjetcms::Real, reco::PattRecoPeak< fftjetcms::Real > > > input_recotree_token_
std::vector< double > doubleBuf
std::vector< fftjet::AbsKernel2d * > memFcns2dVec
std::vector< RecoFFTJet > recoJets
std::unique_ptr< GridAlg > gridAlg
const bool isCrisp
const double recombinationDataCutoff
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleCalcPeak
std::vector< unsigned > occupancy
const double minStableScale
const double maxStableScale
std::vector< fftjet::Peak > preclusters
void beginJob()
Definition: Breakpoints.cc:14
const bool reuseExistingGrid
const unsigned nClustersRequested
std::unique_ptr< std::vector< double > > iniScales
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleCalcJet
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > memberFactorCalcJet
std::unique_ptr< fftjetcms::AbsBgFunctor > bgMembershipFunction
fftjet::AbsVectorRecombinationAlg< fftjetcms::VectorLike, fftjetcms::BgData > RecoAlg
std::unique_ptr< fftjet::Functor2< double, RecoFFTJet, RecoFFTJet > > jetDistanceCalc
math::XYZTLorentzVector VectorLike
unsigned usedLevel
const double noiseLevel
std::vector< RecoFFTJet > iterJets
fftjet::SparseClusteringTree< fftjet::Peak, long > SparseTree
std::vector< unsigned > cellCountsVec
std::unique_ptr< fftjet::ScaleSpaceKernel > jetMembershipFunction
Definition: value.py:1
std::string pileupTableRecord
std::string pileupTableName
std::vector< double > thresholds
fftjet::AbsRecombinationAlg< fftjetcms::Real, fftjetcms::VectorLike, fftjetcms::BgData > GridAlg
const bool calculatePileup
const bool assignConstituents
A grid filled with discretized energy flow.
const edm::InputTag treeLabel
double Real
const edm::InputTag genJetsLabel
const edm::ParameterSet myConfiguration
edm::EDGetTokenT< reco::DiscretizedEnergyFlow > input_energyflow_token_
fftjet::RecombinedJet< fftjetcms::VectorLike > RecoFFTJet
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
const double convergenceDistance
edm::EDGetTokenT< std::vector< reco::FFTAnyJet< reco::GenJet > > > input_genjet_token_
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > recoScaleRatioCalcPeak
const edm::InputTag pileupLabel
Definition: tree.py:1
std::unique_ptr< fftjet::Grid2d< fftjetcms::Real > > pileupEnergyFlow
std::unique_ptr< fftjet::Functor1< double, RecoFFTJet > > recoScaleRatioCalcJet
std::unique_ptr< RecoAlg > recoAlg
fftjetcms::VectorLike unclustered
std::vector< fftjet::Peak > iterPreclusters
std::unique_ptr< fftjet::Functor1< double, fftjet::Peak > > memberFactorCalcPeak
const unsigned maxIterations