CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoLocalCalo/EcalRecProducers/plugins/EcalUncalibRecHitWorkerAnalFit.cc

Go to the documentation of this file.
00001 
00010 #include "RecoLocalCalo/EcalRecProducers/plugins/EcalUncalibRecHitWorkerAnalFit.h"
00011 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00012 #include "DataFormats/EcalDigi/interface/EBDataFrame.h"
00013 #include "DataFormats/EcalDigi/interface/EEDataFrame.h"
00014 #include "DataFormats/EcalDigi/interface/EcalMGPASample.h"
00015 #include "DataFormats/Common/interface/Handle.h"
00016 
00017 #include <iostream>
00018 #include <cmath>
00019 
00020 #include "FWCore/Framework/interface/ESHandle.h"
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 
00024 //#include "CondFormats/EcalObjects/interface/EcalPedestals.h"
00025 //#include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
00026 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00027 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00028 
00029 //#include "CLHEP/Matrix/Matrix.h"
00030 //#include "CLHEP/Matrix/SymMatrix.h"
00031 #include <vector>
00032 
00033 #include "CondFormats/EcalObjects/interface/EcalPedestals.h"
00034 #include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
00035 
00036 #include "CondFormats/EcalObjects/interface/EcalMGPAGainRatio.h"
00037 #include "CondFormats/EcalObjects/interface/EcalGainRatios.h"
00038 #include "CondFormats/DataRecord/interface/EcalGainRatiosRcd.h"
00039 
00040 EcalUncalibRecHitWorkerAnalFit::EcalUncalibRecHitWorkerAnalFit(const edm::ParameterSet& ps) :
00041         EcalUncalibRecHitWorkerBaseClass( ps )
00042 {
00043 }
00044 
00045 
00046 void
00047 EcalUncalibRecHitWorkerAnalFit::set(const edm::EventSetup& es)
00048 {
00049         // Gain Ratios
00050         LogDebug("EcalUncalibRecHitDebug") << "fetching gainRatios....";
00051         es.get<EcalGainRatiosRcd>().get(pRatio);
00052         LogDebug("EcalUncalibRecHitDebug") << "done." ;
00053 
00054         // fetch the pedestals from the cond DB via EventSetup
00055         LogDebug("EcalUncalibRecHitDebug") << "fetching pedestals....";
00056         es.get<EcalPedestalsRcd>().get( pedHandle );
00057         LogDebug("EcalUncalibRecHitDebug") << "done." ;
00058 }
00059 
00060 bool
00061 EcalUncalibRecHitWorkerAnalFit::run( const edm::Event& evt,
00062                 const EcalDigiCollection::const_iterator & itdg,
00063                 EcalUncalibratedRecHitCollection & result)
00064 {
00065         using namespace edm;
00066 
00067         const EcalGainRatioMap & gainMap = pRatio.product()->getMap(); // map of gain ratios
00068         const EcalPedestalsMap & pedMap = pedHandle.product()->getMap(); // map of pedestals
00069 
00070         EcalPedestalsMapIterator pedIter; // pedestal iterator
00071         EcalPedestals::Item aped; // pedestal object for a single xtal
00072 
00073         EcalGainRatioMap::const_iterator gainIter; // gain iterator
00074         EcalMGPAGainRatio aGain; // gain object for a single xtal
00075 
00076         DetId detid( itdg->id() );
00077 
00078         // find pedestals for this channel
00079         //LogDebug("EcalUncalibRecHitDebug") << "looking up pedestal for crystal: " << itdg->id(); // FIXME
00080         pedIter = pedMap.find( detid );
00081         if( pedIter != pedMap.end() ) {
00082                 aped = (*pedIter);
00083         } else {
00084                 edm::LogError("EcalUncalibRecHitWorkerAnalFit") << "error!! could not find pedestals for channel: ";
00085                 if ( detid.subdetId() == EcalBarrel ) {
00086                         edm::LogError("EcalUncalibRecHitWorkerAnalFit") << EBDetId( detid );
00087                 } else {
00088                         edm::LogError("EcalUncalibRecHitWorkerAnalFit") << EEDetId( detid );
00089                 }
00090                 edm::LogError("EcalUncalibRecHitWorkerAnalFit") << "\n  no uncalib rechit will be made for this digi!";
00091                 return false;
00092         }
00093         double pedVec[3];
00094         pedVec[0] = aped.mean_x12;
00095         pedVec[1] = aped.mean_x6;
00096         pedVec[2] = aped.mean_x1;
00097 
00098 
00099         // find gain ratios
00100         //LogDebug("EcalUncalibRecHitDebug") << "looking up gainRatios for crystal: " << itdg->id(); // FIXME
00101         gainIter = gainMap.find( detid );
00102         if( gainIter != gainMap.end() ) {
00103                 aGain = (*gainIter);
00104         } else {
00105                 edm::LogError("EcalUncalibRecHitWorkerAnalFit") << "error!! could not find gain ratios for channel: ";
00106                 if ( detid.subdetId() == EcalBarrel ) {
00107                         edm::LogError("EcalUncalibRecHitWorkerAnalFit") << EBDetId( detid );
00108                 } else {
00109                         edm::LogError("EcalUncalibRecHitWorkerAnalFit") << EEDetId( detid );
00110                 }
00111                 edm::LogError("EcalUncalibRecHitWorkerAnalFit") << "\n  no uncalib rechit will be made for this digi!";
00112                 return false;
00113         }
00114         double gainRatios[3];
00115         gainRatios[0] = 1.;
00116         gainRatios[1] = aGain.gain12Over6();
00117         gainRatios[2] = aGain.gain6Over1()*aGain.gain12Over6();
00118 
00119         if ( detid.subdetId() == EcalBarrel ) {
00120                 EcalUncalibratedRecHit aHit = algoEB_.makeRecHit(*itdg, pedVec, gainRatios, 0 ,0);
00121                 result.push_back( aHit );
00122                 if(aHit.amplitude()>0.) {
00123                         LogDebug("EcalUncalibRecHitInfo") << "EcalUncalibRecHitWorkerAnalFit: processed EBDataFrame with id: "
00124                                 << EBDetId( detid )
00125                                 << "\n" << "uncalib rechit amplitude: " << aHit.amplitude();
00126                 }
00127         } else {
00128                 EcalUncalibratedRecHit aHit = algoEE_.makeRecHit(*itdg, pedVec, gainRatios, 0, 0);
00129                 result.push_back( aHit );
00130                 if(aHit.amplitude()>0.) {
00131                         LogDebug("EcalUncalibRecHitInfo") << "EcalUncalibRecHitWorkerAnalFit: processed EEDataFrame with id: "
00132                                 << EEDetId( detid ) << "\n"
00133                                 << "uncalib rechit amplitude: " << aHit.amplitude();
00134                 }
00135         }
00136         return true;
00137 }
00138 
00139 #include "FWCore/Framework/interface/MakerMacros.h"
00140 #include "RecoLocalCalo/EcalRecProducers/interface/EcalUncalibRecHitWorkerFactory.h"
00141 DEFINE_EDM_PLUGIN( EcalUncalibRecHitWorkerFactory, EcalUncalibRecHitWorkerAnalFit, "EcalUncalibRecHitWorkerAnalFit" );