CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFHBHERecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHBHeRecHitCreator_h
3 
5 
11 
17 
20 public:
22  : PFRecHitCreatorBase(iConfig, cc),
23  recHitToken_(cc.consumes<edm::SortedCollection<HBHERecHit> >(iConfig.getParameter<edm::InputTag>("src"))),
24  geomToken_(cc.esConsumes()) {}
25 
26  void importRecHits(std::unique_ptr<reco::PFRecHitCollection>& out,
27  std::unique_ptr<reco::PFRecHitCollection>& cleaned,
28  const edm::Event& iEvent,
29  const edm::EventSetup& iSetup) override {
30  beginEvent(iEvent, iSetup);
31 
33 
35 
36  // get the ecal geometry
37  const CaloSubdetectorGeometry* hcalBarrelGeo = geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalBarrel);
38  const CaloSubdetectorGeometry* hcalEndcapGeo = geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalEndcap);
39 
40  iEvent.getByToken(recHitToken_, recHitHandle);
41  for (const auto& erh : *recHitHandle) {
42  const HcalDetId detid = erh.idFront();
44 
45  auto energy = erh.energy();
46  auto time = erh.time();
47  auto depth = detid.depth();
48 
49  std::shared_ptr<const CaloCellGeometry> thisCell = nullptr;
51  switch (esd) {
52  case HcalBarrel:
53  thisCell = hcalBarrelGeo->getGeometry(detid);
54  layer = PFLayer::HCAL_BARREL1;
55  break;
56 
57  case HcalEndcap:
58  thisCell = hcalEndcapGeo->getGeometry(detid);
59  layer = PFLayer::HCAL_ENDCAP;
60  break;
61  default:
62  break;
63  }
64 
65  // find rechit geometry
66  if (!thisCell) {
67  edm::LogError("PFHBHERecHitCreator") << "warning detid " << std::hex << detid.rawId() << std::dec << " "
68  << detid << " not found in geometry" << std::endl;
69  continue;
70  }
71 
72  reco::PFRecHit rh(thisCell, detid.rawId(), layer, energy);
73  rh.setTime(time); //Mike: This we will use later
74  rh.setDepth(depth);
75 
76  bool rcleaned = false;
77  bool keep = true;
78 
79  //Apply Q tests
80  for (const auto& qtest : qualityTests_) {
81  if (!qtest->test(rh, erh, rcleaned)) {
82  keep = false;
83  }
84  }
85 
86  if (keep) {
87  out->push_back(std::move(rh));
88  } else if (rcleaned)
89  cleaned->push_back(std::move(rh));
90  }
91  }
92 
93 protected:
95 
96 private:
98 };
99 
100 #endif
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
Log< level::Error, false > LogError
constexpr std::array< uint8_t, layerIndexSize > layer
const int keep
int iEvent
Definition: GenABIO.cc:224
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
def move
Definition: eostools.py:511
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
HcalSubdetector
Definition: HcalAssistant.h:31
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
void setTime(double time)
Definition: PFRecHit.h:73
Layer
layer definition
Definition: PFLayer.h:29
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
PFHBHERecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
void importRecHits(std::unique_ptr< reco::PFRecHitCollection > &out, std::unique_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup) override
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< edm::SortedCollection< HBHERecHit > > recHitToken_