CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
17 
19 
20 template <typename Digi, typename Geometry,PFLayer::Layer Layer,int Detector>
22 
23  public:
25  PFRecHitCreatorBase(iConfig,iC)
26  {
28  }
29 
30  void importRecHits(std::auto_ptr<reco::PFRecHitCollection>&out,std::auto_ptr<reco::PFRecHitCollection>& cleaned ,const edm::Event& iEvent,const edm::EventSetup& iSetup) {
31 
32 
33  beginEvent(iEvent,iSetup);
34 
36 
38  iSetup.get<CaloGeometryRecord>().get(geoHandle);
39 
40  // get the ecal geometry
41  const CaloSubdetectorGeometry *gTmp =
42  geoHandle->getSubdetectorGeometry(DetId::Hcal, Detector);
43 
44  const Geometry *hcalGeo =dynamic_cast< const Geometry* > (gTmp);
45 
46  iEvent.getByToken(recHitToken_,recHitHandle);
47  for( const auto& erh : *recHitHandle ) {
48  const HcalDetId& detid = (HcalDetId)erh.detid();
50 
51  //since hbhe are together kill other detector
52  if (esd !=Detector && Detector != HcalOther )
53  continue;
54 
55 
56  auto energy = erh.energy();
57  auto time = erh.time();
58  auto depth =detid.depth();
59 
60 
61  const CaloCellGeometry * thisCell= hcalGeo->getGeometry(detid);
62 
63  // find rechit geometry
64  if(!thisCell) {
65  edm::LogError("PFHcalRecHitCreator")
66  <<"warning detid "<<detid.rawId()
67  <<" not found in geometry"<<std::endl;
68  continue;
69  }
70 
71 
72  reco::PFRecHit rh(thisCell, detid.rawId(),Layer,
73  energy);
74  rh.setTime(time); //Mike: This we will use later
75  rh.setDepth(depth);
76 
77 
78  bool rcleaned = false;
79  bool keep=true;
80 
81  //Apply Q tests
82  for( const auto& qtest : qualityTests_ ) {
83  if (!qtest->test(rh,erh,rcleaned)) {
84  keep = false;
85 
86  }
87  }
88 
89  if(keep) {
90  out->push_back(std::move(rh));
91  }
92  else if (rcleaned)
93  cleaned->push_back(std::move(rh));
94  }
95  }
96 
97 
98 
99  protected:
101  int hoDepth_;
102 
103 };
104 
110 
111 
112 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< std::unique_ptr< PFRecHitQTestBase > > qualityTests_
PFHcalRecHitCreator< HORecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL2, HcalOuter > PFHORecHitCreator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
void importRecHits(std::auto_ptr< reco::PFRecHitCollection > &out, std::auto_ptr< reco::PFRecHitCollection > &cleaned, const edm::Event &iEvent, const edm::EventSetup &iSetup)
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_BARREL1, HcalBarrel > PFHBRecHitCreator
SeedingLayerSetsHits::SeedingLayer Layer
Definition: LayerTriplets.h:14
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:129
int iEvent
Definition: GenABIO.cc:230
PFHcalRecHitCreator< HBHERecHit, CaloSubdetectorGeometry, PFLayer::HCAL_ENDCAP, HcalEndcap > PFHERecHitCreator
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
def move
Definition: eostools.py:510
HcalSubdetector
Definition: HcalAssistant.h:31
edm::EDGetTokenT< edm::SortedCollection< Digi > > recHitToken_
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
const T & get() const
Definition: EventSetup.h:56
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_EM, HcalForward > PFHFEMRecHitCreator
PFHcalRecHitCreator< HFRecHit, CaloSubdetectorGeometry, PFLayer::HF_HAD, HcalForward > PFHFHADRecHitCreator
PFHcalRecHitCreator(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)