CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCRecHitNavigator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_HGCRecHitNavigator_h
2 #define RecoParticleFlow_PFClusterProducer_HGCRecHitNavigator_h
3 
7 
11 
12 template <ForwardSubdetector D1,
13  typename hgcee,
15  typename hgchef,
17  typename hgcheb>
19 public:
20  HGCRecHitNavigator() = default;
21 
22  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
24 
25  desc.add<std::string>("name", "PFRecHitHGCNavigator");
26 
28  descee.add<std::string>("name", "PFRecHitHGCEENavigator");
29  descee.add<std::string>("topologySource", "HGCalEESensitive");
30  desc.add<edm::ParameterSetDescription>("hgcee", descee);
31 
33  deschef.add<std::string>("name", "PFRecHitHGCHENavigator");
34  deschef.add<std::string>("topologySource", "HGCalHESiliconSensitive");
35  desc.add<edm::ParameterSetDescription>("hgchef", deschef);
36 
38  descheb.add<std::string>("name", "PFRecHitHGCHENavigator");
39  descheb.add<std::string>("topologySource", "HGCalHEScintillatorSensitive");
40  desc.add<edm::ParameterSetDescription>("hgcheb", descheb);
41 
42  descriptions.add("navigator", desc);
43  }
44 
46  if (iConfig.exists("hgcee")) {
47  eeNav_ = new hgcee(iConfig.getParameter<edm::ParameterSet>("hgcee"), cc);
48  } else {
49  eeNav_ = nullptr;
50  }
51  if (iConfig.exists("hgchef")) {
52  hefNav_ = new hgchef(iConfig.getParameter<edm::ParameterSet>("hgchef"), cc);
53  } else {
54  hefNav_ = nullptr;
55  }
56  if (iConfig.exists("hgcheb")) {
57  hebNav_ = new hgcheb(iConfig.getParameter<edm::ParameterSet>("hgcheb"), cc);
58  } else {
59  hebNav_ = nullptr;
60  }
61  }
62 
63  void init(const edm::EventSetup& iSetup) override {
64  if (nullptr != eeNav_)
65  eeNav_->init(iSetup);
66  if (nullptr != hefNav_)
67  hefNav_->init(iSetup);
68  if (nullptr != hebNav_)
69  hebNav_->init(iSetup);
70  }
71 
73  std::unique_ptr<reco::PFRecHitCollection>& hits,
74  edm::RefProd<reco::PFRecHitCollection>& refProd) override {
75  switch (DetId(hit.detId()).subdetId()) {
76  case D1:
77  if (nullptr != eeNav_)
78  eeNav_->associateNeighbours(hit, hits, refProd);
79  break;
80  case D2:
81  if (nullptr != hefNav_)
82  hefNav_->associateNeighbours(hit, hits, refProd);
83  break;
84  case D3:
85  if (nullptr != hebNav_)
86  hebNav_->associateNeighbours(hit, hits, refProd);
87  break;
88  default:
89  break;
90  }
91  }
92 
93 protected:
94  hgcee* eeNav_;
95  hgchef* hefNav_;
96  hgcheb* hebNav_;
97 };
98 
99 #endif
Divides< B, C > D2
Definition: Factorize.h:137
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
bool exists(std::string const &parameterName) const
checks if a parameter exists
HGCRecHitNavigator()=default
ForwardSubdetector
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
void init(const edm::EventSetup &iSetup) override
Divides< A, C > D1
Definition: Factorize.h:136
HGCRecHitNavigator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
void associateNeighbours(reco::PFRecHit &hit, std::unique_ptr< reco::PFRecHitCollection > &hits, edm::RefProd< reco::PFRecHitCollection > &refProd) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: DetId.h:17
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)