00001 #include "RecoLocalCalo/EcalRecProducers/plugins/ESRecHitWorker.h" 00002 00003 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00004 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00005 #include "DataFormats/Common/interface/Handle.h" 00006 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h" 00007 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h" 00008 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h" 00009 #include "CondFormats/DataRecord/interface/ESGainRcd.h" 00010 #include "CondFormats/DataRecord/interface/ESChannelStatusRcd.h" 00011 #include "CondFormats/DataRecord/interface/ESMIPToGeVConstantRcd.h" 00012 #include "CondFormats/DataRecord/interface/ESTimeSampleWeightsRcd.h" 00013 #include "CondFormats/DataRecord/interface/ESPedestalsRcd.h" 00014 #include "CondFormats/DataRecord/interface/ESIntercalibConstantsRcd.h" 00015 #include "CondFormats/DataRecord/interface/ESRecHitRatioCutsRcd.h" 00016 #include "CondFormats/DataRecord/interface/ESAngleCorrectionFactorsRcd.h" 00017 00018 #include <cmath> 00019 #include <iomanip> 00020 #include <iostream> 00021 00022 ESRecHitWorker::ESRecHitWorker(const edm::ParameterSet& ps) : 00023 ESRecHitWorkerBaseClass( ps ) 00024 { 00025 recoAlgo_ = ps.getParameter<int>("ESRecoAlgo"); 00026 00027 if (recoAlgo_ == 0) 00028 algoW_ = new ESRecHitSimAlgo(); 00029 else if (recoAlgo_ == 1) 00030 algoF_ = new ESRecHitFitAlgo(); 00031 else 00032 algoA_ = new ESRecHitAnalyticAlgo(); 00033 } 00034 00035 ESRecHitWorker::~ESRecHitWorker() { 00036 if (recoAlgo_ == 0) 00037 delete algoW_; 00038 else if (recoAlgo_ == 1) 00039 delete algoF_; 00040 else 00041 delete algoA_; 00042 } 00043 00044 void ESRecHitWorker::set(const edm::EventSetup& es) { 00045 00046 es.get<ESGainRcd>().get(esgain_); 00047 const ESGain *gain = esgain_.product(); 00048 00049 es.get<ESMIPToGeVConstantRcd>().get(esMIPToGeV_); 00050 const ESMIPToGeVConstant *mipToGeV = esMIPToGeV_.product(); 00051 00052 double ESGain = gain->getESGain(); 00053 double ESMIPToGeV = (ESGain == 1) ? mipToGeV->getESValueLow() : mipToGeV->getESValueHigh(); 00054 00055 es.get<ESTimeSampleWeightsRcd>().get(esWeights_); 00056 const ESTimeSampleWeights *wgts = esWeights_.product(); 00057 00058 float w0 = wgts->getWeightForTS0(); 00059 float w1 = wgts->getWeightForTS1(); 00060 float w2 = wgts->getWeightForTS2(); 00061 00062 es.get<ESPedestalsRcd>().get(esPedestals_); 00063 const ESPedestals *peds = esPedestals_.product(); 00064 00065 es.get<ESIntercalibConstantsRcd>().get(esMIPs_); 00066 const ESIntercalibConstants *mips = esMIPs_.product(); 00067 00068 es.get<ESAngleCorrectionFactorsRcd>().get(esAngleCorrFactors_); 00069 const ESAngleCorrectionFactors *ang = esAngleCorrFactors_.product(); 00070 00071 es.get<ESChannelStatusRcd>().get(esChannelStatus_); 00072 const ESChannelStatus *channelStatus = esChannelStatus_.product(); 00073 00074 es.get<ESRecHitRatioCutsRcd>().get(esRatioCuts_); 00075 const ESRecHitRatioCuts *ratioCuts = esRatioCuts_.product(); 00076 00077 if (recoAlgo_ == 0) { 00078 algoW_->setESGain(ESGain); 00079 algoW_->setMIPGeV(ESMIPToGeV); 00080 algoW_->setW0(w0); 00081 algoW_->setW1(w1); 00082 algoW_->setW2(w2); 00083 algoW_->setPedestals(peds); 00084 algoW_->setIntercalibConstants(mips); 00085 algoW_->setChannelStatus(channelStatus); 00086 algoW_->setRatioCuts(ratioCuts); 00087 algoW_->setAngleCorrectionFactors(ang); 00088 } else if (recoAlgo_ == 1) { 00089 algoF_->setESGain(ESGain); 00090 algoF_->setMIPGeV(ESMIPToGeV); 00091 algoF_->setPedestals(peds); 00092 algoF_->setIntercalibConstants(mips); 00093 algoF_->setChannelStatus(channelStatus); 00094 algoF_->setRatioCuts(ratioCuts); 00095 algoF_->setAngleCorrectionFactors(ang); 00096 } else { 00097 algoA_->setESGain(ESGain); 00098 algoA_->setMIPGeV(ESMIPToGeV); 00099 algoA_->setPedestals(peds); 00100 algoA_->setIntercalibConstants(mips); 00101 algoA_->setChannelStatus(channelStatus); 00102 algoA_->setRatioCuts(ratioCuts); 00103 algoA_->setAngleCorrectionFactors(ang); 00104 } 00105 } 00106 00107 bool 00108 ESRecHitWorker::run( const ESDigiCollection::const_iterator & itdg, 00109 ESRecHitCollection & result ) 00110 { 00111 if (recoAlgo_ == 0) 00112 result.push_back( algoW_->reconstruct(*itdg) ); 00113 else if (recoAlgo_ == 1) 00114 result.push_back( algoF_->reconstruct(*itdg) ); 00115 else 00116 result.push_back( algoA_->reconstruct(*itdg) ); 00117 return true; 00118 } 00119 00120 #include "FWCore/Framework/interface/MakerMacros.h" 00121 #include "RecoLocalCalo/EcalRecProducers/interface/ESRecHitWorkerFactory.h" 00122 DEFINE_EDM_PLUGIN( ESRecHitWorkerFactory, ESRecHitWorker, "ESRecHitWorker" );