CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFHcalRecHitCreator.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFHcalRecHitCreator_h
2 #define RecoParticleFlow_PFClusterProducer_PFHcalRecHitCreator_h
3 
5 
11 
18 
20 
21 template <typename Digi, typename Geometry, PFLayer::Layer Layer, int Detector>
23 public:
25  : PFRecHitCreatorBase(iConfig, cc),
26  recHitToken_(cc.consumes<edm::SortedCollection<Digi> >(iConfig.getParameter<edm::InputTag>("src"))),
27  geomToken_(cc.esConsumes()),
28  topoToken_(cc.esConsumes()) {}
29 
30  void importRecHits(std::unique_ptr<reco::PFRecHitCollection>& out,
31  std::unique_ptr<reco::PFRecHitCollection>& cleaned,
32  const edm::Event& iEvent,
33  const edm::EventSetup& iSetup) override {
34  beginEvent(iEvent, iSetup);
35 
37 
39  edm::ESHandle<HcalTopology> hcalTopology = iSetup.getHandle(topoToken_);
40 
41  // get the hcal geometry and topology
42  const CaloSubdetectorGeometry* gTmp = geoHandle->getSubdetectorGeometry(DetId::Hcal, Detector);
43  const Geometry* hcalGeo = dynamic_cast<const Geometry*>(gTmp);
44  const HcalTopology* theHcalTopology = hcalTopology.product();
45 
46  iEvent.getByToken(recHitToken_, recHitHandle);
47  for (const auto& erh : *recHitHandle) {
48  HcalDetId detid = (HcalDetId)erh.detid();
50 
51  //since hbhe are together kill other detector
52  if (esd != Detector && Detector != HcalOther)
53  continue;
54 
55  if (theHcalTopology->getMergePositionFlag() && esd == HcalEndcap) {
56  detid = theHcalTopology->idFront(detid);
57  }
58 
59  auto energy = erh.energy();
60  auto time = erh.time();
61  auto depth = detid.depth();
62 
63  auto thisCell = hcalGeo->getGeometry(detid);
64 
65  // find rechit geometry
66  if (!thisCell) {
67  edm::LogError("PFHcalRecHitCreator")
68  << "warning detid " << detid.rawId() << " 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  int hoDepth_;
96 
97 private:
100 };
101 
107 
108 #endif
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
PFHcalRecHitCreator< HORecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL2, HcalOuter > PFHORecHitCreator
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_EM, HcalForward > PFHFEMRecHitCreator
Log< level::Error, false > LogError
bool getMergePositionFlag() const
Definition: HcalTopology.h:167
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
const int keep
int iEvent
Definition: GenABIO.cc:224
Class Geometry Contains vector for fit parameters (mean, sigma, etc.) obtained from multiple IOVs See...
Definition: DMRtrends.cc:183
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
edm::EDGetTokenT< edm::SortedCollection< Digi > > recHitToken_
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_ENDCAP, HcalEndcap > PFHERecHitCreator
void setTime(double time)
Definition: PFRecHit.h:73
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > topoToken_
void importRecHits(std::unique_ptr< reco::PFRecHitCollection > &out, std::unique_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup) override
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_HAD, HcalForward > PFHFHADRecHitCreator
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL1, HcalBarrel > PFHBRecHitCreator
T const * product() const
Definition: ESHandle.h:86
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > geomToken_
PFHcalRecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &cc)
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:170
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283