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  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30 private:
31  // inputs
34  // options
35  // the actual algorithm
36  std::unique_ptr<PFClusterBuilderBase> _pfClusterBuilder;
37  std::unique_ptr<PFClusterEnergyCorrectorBase> _energyCorrector;
38  bool cutsFromDB;
39  HcalPFCuts const* paramPF = nullptr;
40 };
41 
43 
46  desc.add<edm::InputTag>("clustersSource", {});
47  desc.add<edm::ParameterSetDescription>("energyCorrector", {});
48  {
50  pset0.add<std::string>("algoName", "PFMultiDepthClusterizer");
51  {
53  pset1.add<std::string>("algoName", "Basic2DGenericPFlowPositionCalc");
54  {
56  psd.add<std::vector<int>>("depths", {});
57  psd.add<std::string>("detector", "");
58  psd.add<std::vector<double>>("logWeightDenominator", {});
59  pset1.addVPSet("logWeightDenominatorByDetector", psd, {});
60  }
61  pset1.add<double>("minAllowedNormalization", 1e-09);
62  pset1.add<double>("minFractionInCalc", 1e-09);
63  pset1.add<int>("posCalcNCrystals", -1);
64  pset1.add<edm::ParameterSetDescription>("timeResolutionCalcBarrel", {});
65  pset1.add<edm::ParameterSetDescription>("timeResolutionCalcEndcap", {});
66  pset0.add<edm::ParameterSetDescription>("allCellsPositionCalc", pset1);
67  }
68  pset0.add<edm::ParameterSetDescription>("positionCalc", {});
69  pset0.add<double>("minFractionToKeep", 1e-07);
70  pset0.add<double>("nSigmaEta", 2.0);
71  pset0.add<double>("nSigmaPhi", 2.0);
72  desc.add<edm::ParameterSetDescription>("pfClusterBuilder", pset0);
73  }
74  desc.add<edm::ParameterSetDescription>("positionReCalc", {});
75  desc.add<bool>("usePFThresholdsFromDB", false);
76  descriptions.addWithDefaultLabel(desc);
77 }
78 
79 #ifdef PFLOW_DEBUG
80 #define LOGVERB(x) edm::LogVerbatim(x)
81 #define LOGWARN(x) edm::LogWarning(x)
82 #define LOGERR(x) edm::LogError(x)
83 #define LOGDRESSED(x) edm::LogInfo(x)
84 #else
85 #define LOGVERB(x) LogTrace(x)
86 #define LOGWARN(x) edm::LogWarning(x)
87 #define LOGERR(x) edm::LogError(x)
88 #define LOGDRESSED(x) LogDebug(x)
89 #endif
90 
92  _clustersLabel = consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("clustersSource"));
93  cutsFromDB = conf.getParameter<bool>("usePFThresholdsFromDB");
94  const edm::ParameterSet& pfcConf = conf.getParameterSet("pfClusterBuilder");
95 
96  edm::ConsumesCollector&& cc = consumesCollector();
97 
98  if (cutsFromDB) {
99  hcalCutsToken_ = esConsumes<HcalPFCuts, HcalPFCutsRcd, edm::Transition::BeginRun>(edm::ESInputTag("", "withTopo"));
100  }
101 
102  if (!pfcConf.empty()) {
103  const std::string& pfcName = pfcConf.getParameter<std::string>("algoName");
104  if (!pfcName.empty())
105  _pfClusterBuilder = PFClusterBuilderFactory::get()->create(pfcName, pfcConf, cc);
106  }
107  // see if new need to apply corrections, setup if there.
108  const edm::ParameterSet& cConf = conf.getParameterSet("energyCorrector");
109  if (!cConf.empty()) {
110  const std::string& cName = cConf.getParameter<std::string>("algoName");
111  if (!cName.empty())
112  _energyCorrector = PFClusterEnergyCorrectorFactory::get()->create(cName, cConf);
113  }
114 
115  produces<reco::PFClusterCollection>();
116 }
117 
119  if (cutsFromDB) {
121  }
122  _pfClusterBuilder->update(es);
123 }
124 
126  _pfClusterBuilder->reset();
127 
129  e.getByToken(_clustersLabel, inputClusters);
130 
131  std::vector<bool> seedable;
132 
133  auto pfClusters = std::make_unique<reco::PFClusterCollection>();
134  _pfClusterBuilder->buildClusters(*inputClusters, seedable, *pfClusters, paramPF);
135  LOGVERB("PFMultiDepthClusterProducer::produce()") << *_pfClusterBuilder;
136 
137  if (_energyCorrector) {
138  _energyCorrector->correctEnergies(*pfClusters);
139  }
140  e.put(std::move(pfClusters));
141 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
~PFMultiDepthClusterProducer() override=default
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
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
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
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