CMS 3D CMS Logo

HcalHaloDataProducer.cc
Go to the documentation of this file.
3 
4 /*
5  [class]: HcalHaloDataProducer
6  [authors]: R. Remington, The University of Florida
7  [description]: See HcalHaloDataProducer.h
8  [date]: October 15, 2009
9 */
10 
11 using namespace edm;
12 using namespace std;
13 using namespace reco;
14 
15 HcalHaloDataProducer::HcalHaloDataProducer(const edm::ParameterSet& iConfig)
16 {
17  //RecHit Level
18  IT_EBRecHit = iConfig.getParameter<edm::InputTag>("EBRecHitLabel");
19  IT_EERecHit = iConfig.getParameter<edm::InputTag>("EERecHitLabel");
20  IT_HBHERecHit = iConfig.getParameter<edm::InputTag>("HBHERecHitLabel");
21  IT_HFRecHit = iConfig.getParameter<edm::InputTag>("HFRecHitLabel");
22  IT_HORecHit = iConfig.getParameter<edm::InputTag>("HORecHitLabel");
23  IT_CaloTowers = iConfig.getParameter<edm::InputTag>("caloTowerCollName");
24 
25  HBRecHitEnergyThreshold = (float)iConfig.getParameter<double>("HBRecHitEnergyThresholdParam");
26  HERecHitEnergyThreshold = (float)iConfig.getParameter<double>("HERecHitEnergyThresholdParam");
27  SumHcalEnergyThreshold = (float) iConfig.getParameter<double>("SumHcalEnergyThresholdParam");
28  NHitsHcalThreshold = iConfig.getParameter<int>("NHitsHcalThresholdParam");
29 
30  ebrechit_token_ = consumes<EBRecHitCollection>(IT_EBRecHit);
31  eerechit_token_ = consumes<EERecHitCollection>(IT_EERecHit);
32  hbherechit_token_ = consumes<HBHERecHitCollection>(IT_HBHERecHit);
33  hfrechit_token_ = consumes<HFRecHitCollection>(IT_HFRecHit);
34  calotower_token_ = consumes<CaloTowerCollection>(IT_CaloTowers);
35 
36  produces<HcalHaloData>();
37 }
38 
39 void HcalHaloDataProducer::produce(Event& iEvent, const EventSetup& iSetup)
40 {
41  //Get CaloGeometry
42  edm::ESHandle<CaloGeometry> TheCaloGeometry;
43  iSetup.get<CaloGeometryRecord>().get(TheCaloGeometry);
44 
45  //Get CaloTowers
47  iEvent.getByToken(calotower_token_, TheCaloTowers);
48 
49  //Get EB RecHits
51  iEvent.getByToken(ebrechit_token_, TheEBRecHits);
52 
53  //Get EE RecHits
55  iEvent.getByToken(eerechit_token_, TheEERecHits);
56 
57  //Get HB/HE RecHits
58  edm::Handle<HBHERecHitCollection> TheHBHERecHits;
59  // iEvent.getByLabel(IT_HBHERecHit, TheHBHERecHits);
60  iEvent.getByToken(hbherechit_token_, TheHBHERecHits);
61 
62  //Get HF RecHits
64  // iEvent.getByLabel(IT_HFRecHit, TheHFRecHits);
65  iEvent.getByToken(hfrechit_token_, TheHFRecHits);
66 
67  // Run the HcalHaloAlgo to reconstruct the HcalHaloData object
68  HcalHaloAlgo HcalAlgo;
69  HcalAlgo.SetRecHitEnergyThresholds( HBRecHitEnergyThreshold, HERecHitEnergyThreshold );
70  HcalAlgo.SetPhiWedgeThresholds( SumHcalEnergyThreshold, NHitsHcalThreshold );
71 
72 
73 
74  iEvent.put(std::make_unique<HcalHaloData>(HcalAlgo.Calculate(*TheCaloGeometry, TheHBHERecHits, TheCaloTowers, TheEBRecHits, TheEERecHits,iSetup)));
75  return;
76 }
77 
78 HcalHaloDataProducer::~HcalHaloDataProducer(){}
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void SetRecHitEnergyThresholds(float HB, float HE)
Definition: HcalHaloAlgo.h:54
int iEvent
Definition: GenABIO.cc:230
void SetPhiWedgeThresholds(float SumE, int nhits)
Definition: HcalHaloAlgo.h:59
reco::HcalHaloData Calculate(const CaloGeometry &TheCaloGeometry, edm::Handle< HBHERecHitCollection > &TheHBHERecHits, edm::Handle< CaloTowerCollection > &TheCaloTowers, edm::Handle< EBRecHitCollection > &TheEBRecHits, edm::Handle< EERecHitCollection > &TheEERecHits, const edm::EventSetup &TheSetup)
Definition: HcalHaloAlgo.cc:35
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:63