CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFMultiDepthClusterProducer.cc
Go to the documentation of this file.
2 
3 #ifdef PFLOW_DEBUG
4 #define LOGVERB(x) edm::LogVerbatim(x)
5 #define LOGWARN(x) edm::LogWarning(x)
6 #define LOGERR(x) edm::LogError(x)
7 #define LOGDRESSED(x) edm::LogInfo(x)
8 #else
9 #define LOGVERB(x) LogTrace(x)
10 #define LOGWARN(x) edm::LogWarning(x)
11 #define LOGERR(x) edm::LogError(x)
12 #define LOGDRESSED(x) LogDebug(x)
13 #endif
14 
16 {
17  _clustersLabel = consumes<reco::PFClusterCollection>(conf.getParameter<edm::InputTag>("clustersSource"));
18  _pfClusterBuilder.reset(NULL);
19  const edm::ParameterSet& pfcConf = conf.getParameterSet("pfClusterBuilder");
20  if( !pfcConf.empty() ) {
21  const std::string& pfcName = pfcConf.getParameter<std::string>("algoName");
22  PFCBB* pfcb = PFClusterBuilderFactory::get()->create(pfcName,pfcConf);
23  _pfClusterBuilder.reset(pfcb);
24  }
25  // see if new need to apply corrections, setup if there.
26  const edm::ParameterSet& cConf = conf.getParameterSet("energyCorrector");
27  if( !cConf.empty() ) {
28  const std::string& cName = cConf.getParameter<std::string>("algoName");
30  PFClusterEnergyCorrectorFactory::get()->create(cName,cConf);
31  _energyCorrector.reset(eCorr);
32  }
33 
34  produces<reco::PFClusterCollection>();
35 }
36 
38  const edm::EventSetup& es) {
39  _pfClusterBuilder->update(es);
40 
41 }
42 
44  _pfClusterBuilder->reset();
45 
47  e.getByToken(_clustersLabel,inputClusters);
48 
49  std::vector<bool> seedable;
50 
51  std::auto_ptr<reco::PFClusterCollection> pfClusters;
52  pfClusters.reset(new reco::PFClusterCollection);
53  _pfClusterBuilder->buildClusters(*inputClusters, seedable, *pfClusters);
54  LOGVERB("PFMultiDepthClusterProducer::produce()") << *_pfClusterBuilder;
55 
56  if( _energyCorrector ) {
57  _energyCorrector->correctEnergies(*pfClusters);
58  }
59  e.put(pfClusters);
60 
61 }
T getParameter(std::string const &) const
bool empty() const
Definition: ParameterSet.h:218
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
tuple lumi
Definition: fjr2json.py:35
#define NULL
Definition: scimark2.h:8
virtual void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &)
std::unique_ptr< PFClusterBuilderBase > _pfClusterBuilder
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
std::unique_ptr< PFClusterEnergyCorrectorBase > _energyCorrector
tuple conf
Definition: dbtoconf.py:185
PFMultiDepthClusterProducer(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &)
ParameterSet const & getParameterSet(std::string const &) const
edm::EDGetTokenT< reco::PFClusterCollection > _clustersLabel
#define LOGVERB(x)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
static float eCorr(int ieta, int iphi, double ampl, int runnum)
Ugly hack to apply energy corrections to some HB- cells.
T get(const Candidate &c)
Definition: component.h:55