CMS 3D CMS Logo

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 
21  public:
23  PFRecHitCreatorBase(iConfig,iC)
24  {
26  }
27 
28  void importRecHits(std::unique_ptr<reco::PFRecHitCollection>&out,std::unique_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) override {
29 
30  beginEvent(iEvent,iSetup);
31 
33 
35  iSetup.get<CaloGeometryRecord>().get(geoHandle);
36 
37  // get the ecal geometry
38  const CaloSubdetectorGeometry *hcalBarrelGeo = geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalBarrel);
39  const CaloSubdetectorGeometry *hcalEndcapGeo = geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalEndcap);
40 
41  iEvent.getByToken(recHitToken_,recHitHandle);
42  for( const auto& erh : *recHitHandle ) {
43  const HcalDetId detid = erh.idFront();
45 
46  auto energy = erh.energy();
47  auto time = erh.time();
48  auto depth = detid.depth();
49 
50  std::shared_ptr<const CaloCellGeometry> thisCell = nullptr;
52  switch(esd) {
53  case HcalBarrel:
54  thisCell = hcalBarrelGeo->getGeometry(detid);
55  layer = PFLayer::HCAL_BARREL1;
56  break;
57 
58  case HcalEndcap:
59  thisCell = hcalEndcapGeo->getGeometry(detid);
60  layer = PFLayer::HCAL_ENDCAP;
61  break;
62  default:
63  break;
64  }
65 
66  // find rechit geometry
67  if(!thisCell) {
68  edm::LogError("PFHBHERecHitCreator")
69  <<"warning detid "<<std::hex<<detid.rawId()<<std::dec<< " "
70  <<detid<<" not found in geometry"<<std::endl;
71  continue;
72  }
73 
74  reco::PFRecHit rh(thisCell, detid.rawId(),layer,
75  energy);
76  rh.setTime(time); //Mike: This we will use later
77  rh.setDepth(depth);
78 
79  bool rcleaned = false;
80  bool keep=true;
81 
82  //Apply Q tests
83  for( const auto& qtest : qualityTests_ ) {
84  if (!qtest->test(rh,erh,rcleaned)) {
85  keep = false;
86 
87  }
88  }
89 
90  if(keep) {
91  out->push_back(std::move(rh));
92  }
93  else if (rcleaned)
94  cleaned->push_back(std::move(rh));
95  }
96  }
97 
98 
99 
100  protected:
102 
103 
104 };
105 
106 
107 
108 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:44
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const int keep
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
int iEvent
Definition: GenABIO.cc:230
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:32
HcalSubdetector
Definition: HcalAssistant.h:31
void beginEvent(const edm::Event &event, const edm::EventSetup &setup)
void setTime(double time)
Definition: PFRecHit.h:79
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
Layer
layer definition
Definition: PFLayer.h:31
PFHBHERecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
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.
const T & get() const
Definition: EventSetup.h:58
void importRecHits(std::unique_ptr< reco::PFRecHitCollection > &out, std::unique_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup) override
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< edm::SortedCollection< HBHERecHit > > recHitToken_