CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalCalo/EcalRecProducers/plugins/EcalRecHitWorkerSimple.cc

Go to the documentation of this file.
00001 #include "RecoLocalCalo/EcalRecProducers/plugins/EcalRecHitWorkerSimple.h"
00002 
00003 #include "FWCore/Framework/interface/EventSetup.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00006 #include "CondFormats/DataRecord/interface/EcalTimeCalibConstantsRcd.h"
00007 #include "CondFormats/DataRecord/interface/EcalADCToGeVConstantRcd.h"
00008 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
00009 #include "CondFormats/EcalObjects/interface/EcalTimeCalibConstants.h"
00010 #include "CalibCalorimetry/EcalLaserCorrection/interface/EcalLaserDbRecord.h"
00011 
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 
00014 EcalRecHitWorkerSimple::EcalRecHitWorkerSimple(const edm::ParameterSet&ps) :
00015         EcalRecHitWorkerBaseClass(ps)
00016 {
00017         rechitMaker_ = new EcalRecHitSimpleAlgo();
00018         v_chstatus_ = ps.getParameter<std::vector<int> >("ChannelStatusToBeExcluded");
00019         v_DB_reco_flags_ = ps.getParameter<std::vector<int> >("flagsMapDBReco");
00020         killDeadChannels_ = ps.getParameter<bool>("killDeadChannels");
00021         laserCorrection_ = ps.getParameter<bool>("laserCorrection");
00022 }
00023 
00024 
00025 void EcalRecHitWorkerSimple::set(const edm::EventSetup& es)
00026 {
00027         es.get<EcalIntercalibConstantsRcd>().get(ical);
00028         es.get<EcalTimeCalibConstantsRcd>().get(itime);
00029         es.get<EcalADCToGeVConstantRcd>().get(agc);
00030         es.get<EcalChannelStatusRcd>().get(chStatus);
00031         if ( laserCorrection_ ) es.get<EcalLaserDbRecord>().get(laser);
00032 }
00033 
00034 
00035 bool
00036 EcalRecHitWorkerSimple::run( const edm::Event & evt,
00037                 const EcalUncalibratedRecHit& uncalibRH,
00038                 EcalRecHitCollection & result )
00039 {
00040         DetId detid=uncalibRH.id();
00041 
00042         EcalChannelStatusMap::const_iterator chit = chStatus->find(detid);
00043         EcalChannelStatusCode chStatusCode = 1;
00044         if ( chit != chStatus->end() ) {
00045                 chStatusCode = *chit;
00046         } else {
00047                 edm::LogError("EcalRecHitError") << "No channel status found for xtal " 
00048                         << detid.rawId() 
00049                         << "! something wrong with EcalChannelStatus in your DB? ";
00050         }
00051         if ( v_chstatus_.size() > 0) {
00052                 uint16_t code = chStatusCode.getStatusCode() & 0x001F;
00053                 std::vector<int>::const_iterator res = std::find( v_chstatus_.begin(), v_chstatus_.end(), code );
00054                 if ( res != v_chstatus_.end() ) {
00055                         return false;
00056                 }
00057         }
00058 
00059         // find the proper flag for the recHit
00060         // from a configurable vector
00061         // (see cfg file for the association)
00062         uint32_t recoFlag = 0;
00063         uint16_t statusCode = chStatusCode.getStatusCode() & 0x001F;
00064         if ( statusCode < v_DB_reco_flags_.size() ) {
00065                 // not very nice...
00066                 recoFlag = v_DB_reco_flags_[ statusCode ];
00067         } else {
00068                 edm::LogError("EcalRecHitError") << "Flag " << statusCode 
00069                         << " in DB exceed the allowed range of " << v_DB_reco_flags_.size();
00070         }
00071 
00072         const EcalIntercalibConstantMap& icalMap = ical->getMap();  
00073         if ( detid.subdetId() == EcalEndcap ) {
00074                 rechitMaker_->setADCToGeVConstant( float(agc->getEEValue()) );
00075         } else {
00076                 rechitMaker_->setADCToGeVConstant( float(agc->getEBValue()) );
00077         }
00078 
00079         // first intercalibration constants
00080         EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(detid);
00081         EcalIntercalibConstant icalconst = 1;
00082         if( icalit!=icalMap.end() ) {
00083                 icalconst = (*icalit);
00084         } else {
00085                 edm::LogError("EcalRecHitError") << "No intercalib const found for xtal "
00086                         << detid.rawId()
00087                         << "! something wrong with EcalIntercalibConstants in your DB? ";
00088         }
00089 
00090         // get laser coefficient
00091         float lasercalib = 1.;
00092         if ( laserCorrection_ ) lasercalib = laser->getLaserCorrection( detid, evt.time());
00093 
00094         // get time calibration coefficient
00095         const EcalTimeCalibConstantMap & itimeMap = itime->getMap();  
00096         EcalTimeCalibConstantMap::const_iterator itime = itimeMap.find(detid);
00097         EcalTimeCalibConstant itimeconst = 0;
00098         if( itime!=itimeMap.end() ) {
00099                 itimeconst = (*itime);
00100         } else {
00101                 edm::LogError("EcalRecHitError") << "No time calib const found for xtal "
00102                         << detid.rawId()
00103                         << "! something wrong with EcalTimeCalibConstants in your DB? ";
00104         }
00105 
00106         // make the rechit and put in the output collection
00107         if ( recoFlag <= EcalRecHit::kLeadingEdgeRecovered || !killDeadChannels_ ) {
00108                 result.push_back(EcalRecHit( rechitMaker_->makeRecHit(uncalibRH, icalconst * lasercalib, itimeconst, recoFlag) ));
00109         }
00110         return true;
00111 }
00112 
00113 #include "FWCore/Framework/interface/MakerMacros.h"
00114 #include "RecoLocalCalo/EcalRecProducers/interface/EcalRecHitWorkerFactory.h"
00115 DEFINE_EDM_PLUGIN( EcalRecHitWorkerFactory, EcalRecHitWorkerSimple, "EcalRecHitWorkerSimple" );