00001
00002
00003
00004
00005
00021
00022
00023
00024
00025
00026
00027
00028 #ifndef RecoJets_FFTJetProducers_FFTJetProducer_h
00029 #define RecoJets_FFTJetProducers_FFTJetProducer_h
00030
00031 #include <memory>
00032
00033
00034 #include "fftjet/AbsRecombinationAlg.hh"
00035 #include "fftjet/AbsVectorRecombinationAlg.hh"
00036 #include "fftjet/SparseClusteringTree.hh"
00037
00038
00039 #include "FWCore/Framework/interface/Frameworkfwd.h"
00040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00041 #include "FWCore/Framework/interface/EDProducer.h"
00042 #include "FWCore/Framework/interface/Event.h"
00043
00044 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00045
00046
00047 #include "RecoJets/FFTJetAlgorithms/interface/fftjetTypedefs.h"
00048 #include "RecoJets/FFTJetAlgorithms/interface/AbsPileupCalculator.h"
00049 #include "RecoJets/FFTJetProducers/interface/FFTJetInterface.h"
00050
00051 namespace fftjetcms {
00052 class DiscretizedEnergyFlow;
00053 }
00054
00055
00056
00057
00058 class FFTJetProducer : public edm::EDProducer,
00059 public fftjetcms::FFTJetInterface
00060 {
00061 public:
00062 typedef fftjet::RecombinedJet<fftjetcms::VectorLike> RecoFFTJet;
00063 typedef fftjet::SparseClusteringTree<fftjet::Peak,long> SparseTree;
00064
00065
00066 enum StatusBits
00067 {
00068 RESOLUTION = 0xff,
00069 CONSTITUENTS_RESUMMED = 0x100,
00070 PILEUP_CALCULATED = 0x200,
00071 PILEUP_SUBTRACTED_4VEC = 0x400,
00072 PILEUP_SUBTRACTED_PT = 0x800
00073 };
00074
00075 enum Resolution
00076 {
00077 FIXED = 0,
00078 MAXIMALLY_STABLE,
00079 GLOBALLY_ADAPTIVE,
00080 LOCALLY_ADAPTIVE
00081 };
00082
00083 explicit FFTJetProducer(const edm::ParameterSet&);
00084 virtual ~FFTJetProducer();
00085
00086
00087 static Resolution parse_resolution(const std::string& name);
00088
00089 protected:
00090
00091 virtual void beginJob();
00092 virtual void produce(edm::Event&, const edm::EventSetup&);
00093 virtual void endJob();
00094
00095
00096
00097
00098
00099
00100 virtual void selectPreclusters(
00101 const SparseTree& tree,
00102 const fftjet::Functor1<bool,fftjet::Peak>& peakSelector,
00103 std::vector<fftjet::Peak>* preclusters);
00104
00105
00106
00107
00108
00109 virtual void assignMembershipFunctions(
00110 std::vector<fftjet::Peak>* preclusters);
00111
00112
00113 virtual std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> >
00114 parse_peakSelector(const edm::ParameterSet&);
00115
00116
00117 virtual std::auto_ptr<fftjet::ScaleSpaceKernel>
00118 parse_jetMembershipFunction(const edm::ParameterSet&);
00119
00120
00121 virtual std::auto_ptr<fftjetcms::AbsBgFunctor>
00122 parse_bgMembershipFunction(const edm::ParameterSet&);
00123
00124
00125 virtual std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00126 parse_recoScaleCalcPeak(const edm::ParameterSet&);
00127
00128
00129 virtual std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00130 parse_recoScaleRatioCalcPeak(const edm::ParameterSet&);
00131
00132
00133 virtual std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00134 parse_memberFactorCalcPeak(const edm::ParameterSet&);
00135
00136
00137 virtual std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> >
00138 parse_recoScaleCalcJet(const edm::ParameterSet&);
00139
00140 virtual std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> >
00141 parse_recoScaleRatioCalcJet(const edm::ParameterSet&);
00142
00143 virtual std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> >
00144 parse_memberFactorCalcJet(const edm::ParameterSet&);
00145
00146
00147
00148 virtual std::auto_ptr<fftjet::Functor2<double,RecoFFTJet,RecoFFTJet> >
00149 parse_jetDistanceCalc(const edm::ParameterSet&);
00150
00151
00152 virtual std::auto_ptr<fftjetcms::AbsPileupCalculator>
00153 parse_pileupDensityCalc(const edm::ParameterSet& ps);
00154
00155
00156
00157
00158 void selectTreeNodes(const SparseTree& tree,
00159 const fftjet::Functor1<bool,fftjet::Peak>& peakSelect,
00160 std::vector<SparseTree::NodeId>* nodes);
00161 private:
00162 typedef fftjet::AbsVectorRecombinationAlg<
00163 fftjetcms::VectorLike,fftjetcms::BgData> RecoAlg;
00164 typedef fftjet::AbsRecombinationAlg<
00165 fftjetcms::Real,fftjetcms::VectorLike,fftjetcms::BgData> GridAlg;
00166
00167
00168 FFTJetProducer();
00169 FFTJetProducer(const FFTJetProducer&);
00170 FFTJetProducer& operator=(const FFTJetProducer&);
00171
00172
00173 template<class Real>
00174 void loadSparseTreeData(const edm::Event&);
00175
00176 void removeFakePreclusters();
00177
00178
00179
00180 static bool loadEnergyFlow(
00181 const edm::Event& iEvent, const edm::InputTag& label,
00182 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& flow);
00183 void buildGridAlg();
00184 void prepareRecombinationScales();
00185 bool checkConvergence(const std::vector<RecoFFTJet>& previousIterResult,
00186 std::vector<RecoFFTJet>& thisIterResult);
00187 void determineGriddedConstituents();
00188 void determineVectorConstituents();
00189 void saveResults(edm::Event& iEvent, const edm::EventSetup&,
00190 unsigned nPreclustersFound);
00191
00192 template <typename Jet>
00193 void writeJets(edm::Event& iEvent, const edm::EventSetup&);
00194
00195 template <typename Jet>
00196 void makeProduces(const std::string& alias, const std::string& tag);
00197
00198
00199
00200
00201 virtual void determinePileupDensity(
00202 const edm::Event& iEvent, const edm::InputTag& label,
00203 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> >& density);
00204
00205
00206
00207 void determinePileup();
00208
00209
00210
00211
00212
00213
00214
00215 unsigned iterateJetReconstruction();
00216
00217
00218 static void setJetStatusBit(RecoFFTJet* jet, int mask, bool value);
00219
00220
00221
00222
00223
00224
00225 const edm::ParameterSet myConfiguration;
00226
00227
00228 const edm::InputTag treeLabel;
00229
00230
00231
00232 const bool useGriddedAlgorithm;
00233
00234
00235
00236 const bool reuseExistingGrid;
00237
00238
00239 const unsigned maxIterations;
00240
00241
00242 const unsigned nJetsRequiredToConverge;
00243 const double convergenceDistance;
00244
00245
00246 const bool assignConstituents;
00247
00248
00249
00250
00251 const bool resumConstituents;
00252
00253
00254
00255 const bool calculatePileup;
00256 const bool subtractPileup;
00257 const bool subtractPileupAs4Vec;
00258
00259
00260
00261 const edm::InputTag pileupLabel;
00262
00263
00264 const double fixedScale;
00265
00266
00267 const double minStableScale;
00268 const double maxStableScale;
00269
00270
00271 const double stabilityAlpha;
00272
00273
00274
00275 const double noiseLevel;
00276
00277
00278 const unsigned nClustersRequested;
00279
00280
00281 const double gridScanMaxEta;
00282
00283
00284 const std::string recombinationAlgorithm;
00285 const bool isCrisp;
00286 const double unlikelyBgWeight;
00287 const double recombinationDataCutoff;
00288
00289
00290
00291
00292 Resolution resolution;
00293
00294
00295 std::auto_ptr<std::vector<double> > iniScales;
00296
00297
00298 SparseTree sparseTree;
00299
00300
00301 std::auto_ptr<fftjet::Functor1<bool,fftjet::Peak> > peakSelector;
00302
00303
00304 std::auto_ptr<RecoAlg> recoAlg;
00305 std::auto_ptr<GridAlg> gridAlg;
00306 std::auto_ptr<fftjet::ScaleSpaceKernel> jetMembershipFunction;
00307 std::auto_ptr<fftjetcms::AbsBgFunctor> bgMembershipFunction;
00308
00309
00310 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> > recoScaleCalcPeak;
00311
00312
00313 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> >
00314 recoScaleRatioCalcPeak;
00315
00316
00317 std::auto_ptr<fftjet::Functor1<double,fftjet::Peak> > memberFactorCalcPeak;
00318
00319
00320 std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> > recoScaleCalcJet;
00321 std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> > recoScaleRatioCalcJet;
00322 std::auto_ptr<fftjet::Functor1<double,RecoFFTJet> > memberFactorCalcJet;
00323
00324
00325
00326 std::auto_ptr<fftjet::Functor2<double,RecoFFTJet,RecoFFTJet> >
00327 jetDistanceCalc;
00328
00329
00330 std::vector<SparseTree::NodeId> nodes;
00331
00332
00333 std::vector<fftjet::Peak> preclusters;
00334
00335
00336 std::vector<RecoFFTJet> recoJets;
00337
00338
00339 std::vector<unsigned> occupancy;
00340
00341
00342 std::vector<double> thresholds;
00343
00344
00345 unsigned minLevel, maxLevel, usedLevel;
00346
00347
00348 fftjetcms::VectorLike unclustered;
00349 double unused;
00350
00351
00352 std::vector<fftjet::Peak> iterPreclusters;
00353 std::vector<RecoFFTJet> iterJets;
00354 unsigned iterationsPerformed;
00355
00356
00357 std::vector<std::vector<reco::CandidatePtr> > constituents;
00358
00359
00360
00361 std::vector<fftjetcms::VectorLike> pileup;
00362
00363
00364
00365
00366 std::auto_ptr<fftjet::Grid2d<fftjetcms::Real> > pileupEnergyFlow;
00367
00368
00369 std::auto_ptr<fftjetcms::AbsPileupCalculator> pileupDensityCalc;
00370
00371
00372 std::vector<fftjet::AbsKernel2d*> memFcns2dVec;
00373 std::vector<double> doubleBuf;
00374 std::vector<unsigned> cellCountsVec;
00375 };
00376
00377 #endif // RecoJets_FFTJetProducers_FFTJetProducer_h