CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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_HBHERecHit = iConfig.getParameter<edm::InputTag>("HBHERecHitLabel");
19  IT_HFRecHit = iConfig.getParameter<edm::InputTag>("HFRecHitLabel");
20  IT_HORecHit = iConfig.getParameter<edm::InputTag>("HORecHitLabel");
21  IT_CaloTowers = iConfig.getParameter<edm::InputTag>("caloTowerCollName");
22 
23  HBRecHitEnergyThreshold = (float)iConfig.getParameter<double>("HBRecHitEnergyThresholdParam");
24  HERecHitEnergyThreshold = (float)iConfig.getParameter<double>("HERecHitEnergyThresholdParam");
25  SumHcalEnergyThreshold = (float) iConfig.getParameter<double>("SumHcalEnergyThresholdParam");
26  NHitsHcalThreshold = iConfig.getParameter<int>("NHitsHcalThresholdParam");
27 
28  hbherechit_token_ = consumes<HBHERecHitCollection>(IT_HBHERecHit);
29  hfrechit_token_ = consumes<HFRecHitCollection>(IT_HFRecHit);
30  calotower_token_ = consumes<CaloTowerCollection>(IT_CaloTowers);
31 
32  produces<HcalHaloData>();
33 }
34 
35 void HcalHaloDataProducer::produce(Event& iEvent, const EventSetup& iSetup)
36 {
37  //Get CaloGeometry
38  edm::ESHandle<CaloGeometry> TheCaloGeometry;
39  iSetup.get<CaloGeometryRecord>().get(TheCaloGeometry);
40 
41  //Get CaloTowers
43  iEvent.getByToken(calotower_token_, TheCaloTowers);
44 
45  //Get HB/HE RecHits
46  edm::Handle<HBHERecHitCollection> TheHBHERecHits;
47  // iEvent.getByLabel(IT_HBHERecHit, TheHBHERecHits);
48  iEvent.getByToken(hbherechit_token_, TheHBHERecHits);
49 
50  //Get HF RecHits
52  // iEvent.getByLabel(IT_HFRecHit, TheHFRecHits);
53  iEvent.getByToken(hfrechit_token_, TheHFRecHits);
54 
55  // Run the HcalHaloAlgo to reconstruct the HcalHaloData object
56  HcalHaloAlgo HcalAlgo;
57  HcalAlgo.SetRecHitEnergyThresholds( HBRecHitEnergyThreshold, HERecHitEnergyThreshold );
58  HcalAlgo.SetPhiWedgeThresholds( SumHcalEnergyThreshold, NHitsHcalThreshold );
59 
60  HcalHaloData HcalData;
61  if( TheCaloGeometry.isValid() && TheHBHERecHits.isValid() )
62  {
63  if( TheCaloTowers.isValid() ) {
64  std::auto_ptr<HcalHaloData> HcalData( new HcalHaloData( HcalAlgo.Calculate(*TheCaloGeometry, TheHBHERecHits, TheCaloTowers) ) ) ;
65  iEvent.put ( HcalData ) ;
66  } else {
67  std::auto_ptr<HcalHaloData> HcalData( new HcalHaloData( HcalAlgo.Calculate(*TheCaloGeometry, TheHBHERecHits) ) ) ;
68  iEvent.put ( HcalData ) ;
69  }
70  }
71  else
72  {
73  std::auto_ptr<HcalHaloData> HcalData( new HcalHaloData() ) ;
74  iEvent.put( HcalData ) ;
75  }
76  return;
77 }
78 
79 HcalHaloDataProducer::~HcalHaloDataProducer(){}
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void SetRecHitEnergyThresholds(float HB, float HE)
Definition: HcalHaloAlgo.h:43
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
reco::HcalHaloData Calculate(const CaloGeometry &TheCaloGeometry, edm::Handle< HBHERecHitCollection > &TheHBHERecHits, edm::Handle< CaloTowerCollection > &TheCaloTowers)
Definition: HcalHaloAlgo.cc:33
bool isValid() const
Definition: HandleBase.h:75
void SetPhiWedgeThresholds(float SumE, int nhits)
Definition: HcalHaloAlgo.h:48
const T & get() const
Definition: EventSetup.h:56
bool isValid() const
Definition: ESHandle.h:47