#include <EcalUncalibRecHitWorkerFixedAlphaBetaFit.h>
Public Member Functions | |
EcalUncalibRecHitWorkerFixedAlphaBetaFit (const edm::ParameterSet &ps) | |
bool | run (const edm::Event &evt, const EcalDigiCollection::const_iterator &digi, EcalUncalibratedRecHitCollection &result) |
void | set (const edm::EventSetup &es) |
virtual | ~EcalUncalibRecHitWorkerFixedAlphaBetaFit () |
Private Member Functions | |
bool | setAlphaBeta () |
Private Attributes | |
EcalUncalibRecHitFixedAlphaBetaAlgo < EBDataFrame > | algoEB_ |
EcalUncalibRecHitFixedAlphaBetaAlgo < EEDataFrame > | algoEE_ |
std::string | alphabetaFilename_ |
std::vector< std::vector < std::pair< double, double > > > | alphaBetaValues_ |
double | alphaEB_ |
double | alphaEE_ |
double | AmplThrEB_ |
double | AmplThrEE_ |
double | betaEB_ |
double | betaEE_ |
edm::ESHandle< EcalPedestals > | pedHandle |
edm::ESHandle< EcalGainRatios > | pRatio |
bool | useAlphaBetaArray_ |
Definition at line 22 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
EcalUncalibRecHitWorkerFixedAlphaBetaFit::EcalUncalibRecHitWorkerFixedAlphaBetaFit | ( | const edm::ParameterSet & | ps | ) |
Definition at line 31 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.cc.
References algoEB_, algoEE_, alphabetaFilename_, alphaEB_, alphaEE_, betaEB_, betaEE_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), setAlphaBeta(), EcalUncalibRecHitFixedAlphaBetaAlgo< C >::SetDynamicPedestal(), EcalUncalibRecHitFixedAlphaBetaAlgo< C >::SetMinAmpl(), and useAlphaBetaArray_.
: EcalUncalibRecHitWorkerBaseClass( ps ) { alphaEB_= ps.getParameter<double>("alphaEB"); betaEB_= ps.getParameter<double>("betaEB"); alphaEE_= ps.getParameter<double>("alphaEE"); betaEE_= ps.getParameter<double>("betaEE"); alphabetaFilename_= ps.getUntrackedParameter<std::string>("AlphaBetaFilename","NOFILE"); useAlphaBetaArray_=setAlphaBeta(); // set crystalwise values of alpha and beta if ( !useAlphaBetaArray_ ) { edm::LogInfo("EcalUncalibRecHitError") << " No alfa-beta file found. Using the deafult values."; } algoEB_.SetMinAmpl( ps.getParameter<double> ("MinAmplBarrel") ); algoEE_.SetMinAmpl( ps.getParameter<double> ("MinAmplEndcap") ); bool dyn_pede = ps.getParameter<bool>("UseDynamicPedestal"); algoEB_.SetDynamicPedestal(dyn_pede); algoEE_.SetDynamicPedestal(dyn_pede); }
virtual EcalUncalibRecHitWorkerFixedAlphaBetaFit::~EcalUncalibRecHitWorkerFixedAlphaBetaFit | ( | ) | [inline, virtual] |
Definition at line 26 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
{};
bool EcalUncalibRecHitWorkerFixedAlphaBetaFit::run | ( | const edm::Event & | evt, |
const EcalDigiCollection::const_iterator & | digi, | ||
EcalUncalibratedRecHitCollection & | result | ||
) | [virtual] |
Implements EcalUncalibRecHitWorkerBaseClass.
Definition at line 104 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.cc.
References a, algoEB_, algoEE_, alphaBetaValues_, alphaEB_, alphaEE_, b, betaEB_, betaEE_, cond::rpcobgas::detid, EcalBarrel, EcalCondObjectContainer< T >::end(), EcalCondObjectContainer< T >::find(), EcalMGPAGainRatio::gain12Over6(), EcalMGPAGainRatio::gain6Over1(), EBDetId::ic(), EBDetId::ism(), EcalUncalibRecHitFixedAlphaBetaAlgo< C >::makeRecHit(), pedHandle, pRatio, edm::ESHandle< T >::product(), edm::SortedCollection< T, SORT >::push_back(), EcalUncalibRecHitFixedAlphaBetaAlgo< C >::SetAlphaBeta(), and useAlphaBetaArray_.
{ const EcalGainRatioMap& gainMap = pRatio.product()->getMap(); // map of gain ratios EcalGainRatioMap::const_iterator gainIter; // gain iterator EcalMGPAGainRatio aGain; // gain object for a single xtal const EcalPedestalsMap & pedMap = pedHandle.product()->getMap(); // map of pedestals EcalPedestalsMapIterator pedIter; // pedestal iterator EcalPedestals::Item aped; // pedestal object for a single xtal DetId detid( itdg->id() ); // find pedestals for this channel //LogDebug("EcalUncalibRecHitDebug") << "looking up pedestal for crystal: " << itdg->id(); pedIter = pedMap.find(itdg->id()); if( pedIter != pedMap.end() ) { aped = (*pedIter); } else { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << "error!! could not find pedestals for channel: "; if ( detid.subdetId() == EcalBarrel ) { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << EBDetId( detid ); } else { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << EEDetId( detid ); } edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << "\n no uncalib rechit will be made for this digi!"; return false; } double pedVec[3]; pedVec[0] = aped.mean_x12; pedVec[1] = aped.mean_x6; pedVec[2] = aped.mean_x1; // find gain ratios //LogDebug("EcalUncalibRecHitDebug") << "looking up gainRatios for crystal: " << EBDetId(itdg->id()) ; // FIXME!!!!!!!! gainIter = gainMap.find(itdg->id()); if( gainIter != gainMap.end() ) { aGain = (*gainIter); } else { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << "error!! could not find gain ratios for channel: "; if ( detid.subdetId() == EcalBarrel ) { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << EBDetId( detid ); } else { edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << EEDetId( detid ); } edm::LogError("EcalUncalibRecHitWorkerFixedAlphaBetaFit") << "\n no uncalib rechit will be made for this digi!"; return false; } double gainRatios[3]; gainRatios[0] = 1.; gainRatios[1] = aGain.gain12Over6(); gainRatios[2] = aGain.gain6Over1()*aGain.gain12Over6(); if ( detid.subdetId() == EcalBarrel ) { // Define Alpha and Beta either by stored values or by default universal values EBDetId ebDetId( detid ); double a, b; if (useAlphaBetaArray_){ if ( alphaBetaValues_[ ebDetId.ism()-1 ].size() != 0 ) { a = alphaBetaValues_[ebDetId.ism()-1][ebDetId.ic()-1].first; b = alphaBetaValues_[ebDetId.ism()-1][ebDetId.ic()-1].second; if ( ( a == -1 ) && ( b == -1 ) ) { a = alphaEB_; b = betaEB_; } } else { a = alphaEB_; b = betaEB_; } } else { a = alphaEB_; b = betaEB_; } algoEB_.SetAlphaBeta(a,b); result.push_back( algoEB_.makeRecHit( *itdg, pedVec, gainRatios, 0, 0) ); } else { //FIX ME load in a and b from a file algoEE_.SetAlphaBeta(alphaEE_,betaEE_); result.push_back( algoEE_.makeRecHit(*itdg, pedVec, gainRatios, 0 , 0) ); } return true; }
void EcalUncalibRecHitWorkerFixedAlphaBetaFit::set | ( | const edm::EventSetup & | es | ) | [virtual] |
Implements EcalUncalibRecHitWorkerBaseClass.
Definition at line 55 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.cc.
References edm::EventSetup::get(), LogDebug, pedHandle, and pRatio.
{ // Gain Ratios LogDebug("EcalUncalibRecHitDebug") << "fetching gainRatios...."; es.get<EcalGainRatiosRcd>().get(pRatio); LogDebug("EcalUncalibRecHitDebug") << "done." ; // fetch the pedestals from the cond DB via EventSetup LogDebug("EcalUncalibRecHitDebug") << "fetching pedestals...."; es.get<EcalPedestalsRcd>().get( pedHandle ); LogDebug("EcalUncalibRecHitDebug") << "done." ; }
bool EcalUncalibRecHitWorkerFixedAlphaBetaFit::setAlphaBeta | ( | ) | [private] |
Definition at line 70 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.cc.
References a, alphabetaFilename_, alphaBetaValues_, b, mergeVDriftHistosByStation::file, AlCaHLTBitMon_ParallelJobs::p, runTheMatrix::ret, and findQualityFiles::size.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit().
{ std::ifstream file(alphabetaFilename_.c_str()); if (! file.is_open()) return false; alphaBetaValues_.resize(36); char buffer[100]; int sm, cry,ret; float a,b; std::pair<double,double> p(-1,-1); while( ! file.getline(buffer,100).eof() ){ ret=sscanf(buffer,"%d %d %f %f", &sm, &cry, &a, &b); if ((ret!=4)|| (sm<=0) ||(sm>36)|| (cry<=0)||(cry>1700)){ // send warning continue; } if (alphaBetaValues_[sm-1].size()==0){ alphaBetaValues_[sm-1].resize(1700,p); } alphaBetaValues_[sm-1][cry-1].first = a; alphaBetaValues_[sm-1][cry-1].second = b; } file.close(); return true; }
EcalUncalibRecHitFixedAlphaBetaAlgo<EBDataFrame> EcalUncalibRecHitWorkerFixedAlphaBetaFit::algoEB_ [private] |
Definition at line 36 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
EcalUncalibRecHitFixedAlphaBetaAlgo<EEDataFrame> EcalUncalibRecHitWorkerFixedAlphaBetaFit::algoEE_ [private] |
Definition at line 37 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
std::string EcalUncalibRecHitWorkerFixedAlphaBetaFit::alphabetaFilename_ [private] |
Definition at line 45 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and setAlphaBeta().
std::vector<std::vector<std::pair<double,double> > > EcalUncalibRecHitWorkerFixedAlphaBetaFit::alphaBetaValues_ [private] |
Definition at line 43 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by run(), and setAlphaBeta().
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::alphaEB_ [private] |
Definition at line 39 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::alphaEE_ [private] |
Definition at line 41 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::AmplThrEB_ [private] |
Definition at line 33 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::AmplThrEE_ [private] |
Definition at line 34 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::betaEB_ [private] |
Definition at line 40 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
double EcalUncalibRecHitWorkerFixedAlphaBetaFit::betaEE_ [private] |
Definition at line 42 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().
Definition at line 50 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Definition at line 49 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
bool EcalUncalibRecHitWorkerFixedAlphaBetaFit::useAlphaBetaArray_ [private] |
Definition at line 44 of file EcalUncalibRecHitWorkerFixedAlphaBetaFit.h.
Referenced by EcalUncalibRecHitWorkerFixedAlphaBetaFit(), and run().