CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoLocalCalo/HcalRecProducers/src/HcalSimpleReconstructor.cc

Go to the documentation of this file.
00001 #include "HcalSimpleReconstructor.h"
00002 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00003 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00004 #include "DataFormats/Common/interface/EDCollection.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/Selector.h"
00007 #include "FWCore/Framework/interface/ESHandle.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00010 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00011 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00012 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00013 
00014 #include <iostream>
00015     
00016 HcalSimpleReconstructor::HcalSimpleReconstructor(edm::ParameterSet const& conf):
00017   reco_(conf.getParameter<bool>("correctForTimeslew"),
00018         conf.getParameter<bool>("correctForPhaseContainment"),conf.getParameter<double>("correctionPhaseNS")),
00019   det_(DetId::Hcal),
00020   inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
00021   dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
00022   firstSample_(conf.getParameter<int>("firstSample")),
00023   samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
00024   tsFromDB_(conf.getParameter<bool>("tsFromDB"))
00025 {
00026 
00027   std::string subd=conf.getParameter<std::string>("Subdetector");
00028   if (!strcasecmp(subd.c_str(),"HBHE")) {
00029     subdet_=HcalBarrel;
00030     produces<HBHERecHitCollection>();
00031   } else if (!strcasecmp(subd.c_str(),"HO")) {
00032     subdet_=HcalOuter;
00033     produces<HORecHitCollection>();
00034   } else if (!strcasecmp(subd.c_str(),"HF")) {
00035     subdet_=HcalForward;
00036     produces<HFRecHitCollection>();
00037   } else {
00038     std::cout << "HcalSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
00039   }       
00040   
00041 }
00042 
00043 HcalSimpleReconstructor::~HcalSimpleReconstructor() { }
00044 
00045 void HcalSimpleReconstructor::beginRun(edm::Run&r, edm::EventSetup const & es){
00046   if(tsFromDB_) {
00047     edm::ESHandle<HcalRecoParams> p;
00048     es.get<HcalRecoParamsRcd>().get(p);
00049     paramTS = new HcalRecoParams(*p.product());
00050   }
00051   reco_.beginRun(es);
00052 }
00053 
00054 void HcalSimpleReconstructor::endRun(edm::Run&r, edm::EventSetup const & es){
00055   if(tsFromDB_ && paramTS) {
00056     delete paramTS;
00057     paramTS = 0;
00058     reco_.endRun();
00059   }
00060 }
00061 
00062 
00063 template<class DIGICOLL, class RECHITCOLL> 
00064 void HcalSimpleReconstructor::process(edm::Event& e, const edm::EventSetup& eventSetup)
00065 {
00066   // get conditions
00067   edm::ESHandle<HcalDbService> conditions;
00068   eventSetup.get<HcalDbRecord>().get(conditions);
00069   const HcalQIEShape* shape = conditions->getHcalShape (); // this one is generic
00070 
00071   // HACK related to HB- corrections
00072   if(e.isRealData()) reco_.setForData();
00073 
00074   edm::Handle<DIGICOLL> digi;
00075 
00076   e.getByLabel(inputLabel_,digi);
00077 
00078   // create empty output
00079   std::auto_ptr<RECHITCOLL> rec(new RECHITCOLL);
00080   rec->reserve(digi->size());
00081   // run the algorithm
00082   int first = firstSample_;
00083   int toadd = samplesToAdd_;
00084   typename DIGICOLL::const_iterator i;
00085   for (i=digi->begin(); i!=digi->end(); i++) {
00086     HcalDetId cell = i->id();
00087     DetId detcell=(DetId)cell;
00088     // rof 27.03.09: drop ZS marked and passed digis:
00089     if (dropZSmarkedPassed_)
00090       if (i->zsMarkAndPass()) continue;
00091 
00092     const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
00093     const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00094     HcalCoderDb coder (*channelCoder, *shape);
00095 
00096     //>>> firstSample & samplesToAdd
00097     if(tsFromDB_) {
00098       const HcalRecoParam* param_ts = paramTS->getValues(detcell.rawId());
00099       first = param_ts->firstSample();
00100       toadd = param_ts->samplesToAdd();
00101     }
00102     rec->push_back(reco_.reconstruct(*i,first,toadd,coder,calibrations));
00103 
00104   }
00105   // return result
00106   e.put(rec);
00107 }
00108 
00109 
00110 void HcalSimpleReconstructor::produce(edm::Event& e, const edm::EventSetup& eventSetup)
00111 {
00112   // HACK related to HB- corrections
00113   if(e.isRealData()) reco_.setForData();
00114  
00115   
00116   if (det_==DetId::Hcal) {
00117     if (subdet_==HcalBarrel || subdet_==HcalEndcap) {
00118       process<HBHEDigiCollection, HBHERecHitCollection>(e, eventSetup);
00119     } else if (subdet_==HcalForward) {
00120       process<HFDigiCollection, HFRecHitCollection>(e, eventSetup);
00121     } else if (subdet_==HcalOuter) {
00122       process<HODigiCollection, HORecHitCollection>(e, eventSetup);
00123     } else if (subdet_==HcalOther && subdetOther_==HcalCalibration) {
00124       process<HcalCalibDigiCollection, HcalCalibRecHitCollection>(e, eventSetup);
00125     }
00126   } 
00127 }