CMS 3D CMS Logo

PFMultiDepthClusterProducer.cc
Go to the documentation of this file.
14 
15 #include <memory>
16 
21 
22 public:
24  ~PFMultiDepthClusterProducer() override = default;
25 
26  void beginRun(const edm::Run&, const edm::EventSetup&) override;
27  void produce(edm::Event&, const edm::EventSetup&) override;
28 
29 private:
30  // inputs
33  // options
34  // the actual algorithm
35  std::unique_ptr<PFClusterBuilderBase> _pfClusterBuilder;
36  std::unique_ptr<PFClusterEnergyCorrectorBase> _energyCorrector;
37  bool cutsFromDB;
38  HcalPFCuts const* paramPF = nullptr;
39 };
40 
42 
43 #ifdef PFLOW_DEBUG
44 #define LOGVERB(x) edm::LogVerbatim(x)
45 #define LOGWARN(x) edm::LogWarning(x)
46 #define LOGERR(x) edm::LogError(x)
47 #define LOGDRESSED(x) edm::LogInfo(x)
48 #else
49 #define LOGVERB(x) LogTrace(x)
50 #define LOGWARN(x) edm::LogWarning(x)
51 #define LOGERR(x) edm::LogError(x)
52 #define LOGDRESSED(x) LogDebug(x)
53 #endif
54 
56  _clustersLabel = consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("clustersSource"));
57  cutsFromDB = conf.getParameter<bool>("usePFThresholdsFromDB");
58  const edm::ParameterSet& pfcConf = conf.getParameterSet("pfClusterBuilder");
59 
60  edm::ConsumesCollector&& cc = consumesCollector();
61 
62  if (cutsFromDB) {
63  hcalCutsToken_ = esConsumes<HcalPFCuts, HcalPFCutsRcd, edm::Transition::BeginRun>(edm::ESInputTag("", "withTopo"));
64  }
65 
66  if (!pfcConf.empty()) {
67  const std::string& pfcName = pfcConf.getParameter<std::string>("algoName");
68  _pfClusterBuilder = PFClusterBuilderFactory::get()->create(pfcName, pfcConf, cc);
69  }
70  // see if new need to apply corrections, setup if there.
71  const edm::ParameterSet& cConf = conf.getParameterSet("energyCorrector");
72  if (!cConf.empty()) {
73  const std::string& cName = cConf.getParameter<std::string>("algoName");
75  }
76 
77  produces<reco::PFClusterCollection>();
78 }
79 
81  if (cutsFromDB) {
83  }
84  _pfClusterBuilder->update(es);
85 }
86 
88  _pfClusterBuilder->reset();
89 
91  e.getByToken(_clustersLabel, inputClusters);
92 
93  std::vector<bool> seedable;
94 
95  auto pfClusters = std::make_unique<reco::PFClusterCollection>();
96  _pfClusterBuilder->buildClusters(*inputClusters, seedable, *pfClusters, paramPF);
97  LOGVERB("PFMultiDepthClusterProducer::produce()") << *_pfClusterBuilder;
98 
99  if (_energyCorrector) {
100  _energyCorrector->correctEnergies(*pfClusters);
101  }
102  e.put(std::move(pfClusters));
103 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
~PFMultiDepthClusterProducer() override=default
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
ParameterSet const & getParameterSet(std::string const &) const
std::unique_ptr< PFClusterBuilderBase > _pfClusterBuilder
bool empty() const
Definition: ParameterSet.h:202
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::unique_ptr< PFClusterEnergyCorrectorBase > _energyCorrector
edm::ESGetToken< HcalPFCuts, HcalPFCutsRcd > hcalCutsToken_
void produce(edm::Event &, const edm::EventSetup &) override
PFMultiDepthClusterProducer(const edm::ParameterSet &)
edm::EDGetTokenT< reco::PFClusterCollection > _clustersLabel
#define LOGVERB(x)
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define get
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45