HGCal Test Beam  03a93d6239a951948e06fb3ef8dae4cbdebfad30
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
HGCalTBRecHitProducer.cc
Go to the documentation of this file.
1 #include "HGCal/Reco/plugins/HGCalTBRecHitProducer.h"
2 
3 HGCalTBRecHitProducer::HGCalTBRecHitProducer(const edm::ParameterSet& cfg)
4  : outputCollectionName(cfg.getParameter<std::string>("OutputCollectionName")),
5  _digisToken(consumes<HGCalTBDigiCollection>(cfg.getParameter<edm::InputTag>("digiCollection"))),
6  _pedestalLow_filename(cfg.getParameter<std::string>("pedestalLow")),
7  _pedestalHigh_filename(cfg.getParameter<std::string>("pedestalHigh")),
8  _gainsLow_filename(cfg.getParameter<std::string>("gainLow")),
9  _gainsHigh_filename(cfg.getParameter<std::string>("gainHigh")),
10  _adcSaturation(cfg.getParameter<unsigned int>("adcSaturation"))
11 {
12  produces <HGCalTBRecHitCollection>(outputCollectionName);
13 
14 }
15 
16 void HGCalTBRecHitProducer::produce(edm::Event& event, const edm::EventSetup& iSetup)
17 {
18 
19  std::auto_ptr<HGCalTBRecHitCollection> rechits(new HGCalTBRecHitCollection); //auto_ptr are automatically deleted when out of scope
20 
21  edm::Handle<HGCalTBDigiCollection> digisHandle;
22  event.getByToken(_digisToken, digisHandle);
23 
24 #ifdef ESPRODUCER
25 // the conditions should be defined by ESProduers and available in the iSetup
26  edm::ESHandle<HGCalCondPedestals> pedestalsHandle;
27  //iSetup.get<HGCalCondPedestals>().get(pedestalsHandle); ///\todo should we support such a syntax?
28  HGCalCondPedestals pedestals = iSetup.get<HGCalCondPedestals>();
29 
30  edm::ESHandle<HGCalCondGains> adcToGeVHandle;
31  //iSetup.get<HGCalCondGains>().get(adcToGeVHandle);
32  HGCalCondGains adcToGeV = iSetup.get<HGCalCondGains>();
33 #else
34  HGCalCondGains adcToGeV_low, adcToGeV_high;
35  HGCalCondPedestals pedestals_low, pedestals_high;
37  //HGCalElectronicsMap emap;
38 // assert(io.load("mapfile.txt",emap)); ///\todo to be trasformed into exception
39  if(_pedestalLow_filename != "") assert(condIO.load(_pedestalLow_filename, pedestals_low));
40  if(_pedestalHigh_filename != "") assert(condIO.load(_pedestalHigh_filename, pedestals_high));
41  if(_gainsLow_filename != "") assert(condIO.load(_gainsLow_filename, adcToGeV_low));
42  if(_gainsHigh_filename != "") assert(condIO.load(_gainsHigh_filename, adcToGeV_high));
43  ///\todo check if reading the conditions from file some channels are not in the file!
44 #endif
45 
46  for(auto digi_itr = digisHandle->begin(); digi_itr != digisHandle->end(); ++digi_itr) {
47 #ifdef DEBUG
48  std::cout << "[RECHIT PRODUCER: digi]" << *digi_itr << std::endl;
49 #endif
50 
51 //------------------------------ this part should go into HGCalTBRecoAlgo class
52 
53  SKIROC2DataFrame digi(*digi_itr);
54  unsigned int nSamples = digi.samples();
55 
56  // if there are more than 1 sample, we need to define a reconstruction algorithm
57  // now taking the first sample
58  for(unsigned int iSample = 0; iSample < nSamples; ++iSample) {
59 
60  float pedestal_low_value = (pedestals_low.size() > 0) ? pedestals_low.get(digi.detid())->value : 0;
61  float pedestal_high_value = (pedestals_high.size() > 0) ? pedestals_high.get(digi.detid())->value : 0;
62  float adcToGeV_low_value = (adcToGeV_low.size() > 0) ? adcToGeV_low.get(digi.detid())->value : 1;
63  float adcToGeV_high_value = (adcToGeV_high.size() > 0) ? adcToGeV_high.get(digi.detid())->value : 1;
64 
65  float energyLow = digi[iSample].adcLow() - pedestal_low_value * adcToGeV_low_value;
66  float energyHigh = digi[iSample].adcHigh() - pedestal_high_value * adcToGeV_high_value;
67 
68  HGCalTBRecHit recHit(digi.detid(), energyLow, energyHigh, digi[iSample].tdc()); ///\todo use time calibration!
69 
70  if(digi[iSample].adcHigh() > _adcSaturation) recHit.setFlag(kHighGainSaturated);
71  if(digi[iSample].adcLow() > _adcSaturation) recHit.setFlag(kLowGainSaturated);
72 
73 
74 #ifdef DEBUG
75  std::cout << recHit << std::endl;
76 #endif
77  if(iSample == 0) rechits->push_back(recHit); ///\todo define an algorithm for the energy if more than 1 sample, code inefficient
78  }
79  }
80  event.put(rechits, outputCollectionName);
81 }
82 // Should there be a destructor ??
kLowGainSaturated
Definition: HGCalTBRecHit.h:20
DEFINE_FWK_MODULE(Pedestals)
void setFlag(int flag)
Definition: HGCalTBRecHit.h:50
kHighGainSaturated
Definition: HGCalTBRecHit.h:20
edm::SortedCollection< HGCalTBRecHit > HGCalTBRecHitCollection
HGCalTBDetId detid() const
Get the detector id.
static const HGCalCondObjectNumberingScheme * scheme()
bool load(const std::string &filename, HGCalCondObjectContainer< float > &)
load conditions from file
const Item * get(DetId id) const
HGCalTBRecHitProducer(const edm::ParameterSet &)
virtual void produce(edm::Event &, const edm::EventSetup &)
int samples() const
total number of samples in the digi