CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/CalibFormats/HcalObjects/src/HcalDbService.cc

Go to the documentation of this file.
00001 //
00002 // F.Ratnikov (UMd), Aug. 9, 2005
00003 //
00004 
00005 #include "FWCore/Utilities/interface/typelookup.h"
00006 
00007 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00008 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00009 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00010 #include "CalibFormats/HcalObjects/interface/HcalCalibrationWidths.h"
00011 
00012 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
00013 
00014 #include <cmath>
00015 
00016 HcalDbService::HcalDbService (const edm::ParameterSet& cfg): 
00017   mPedestals (0), mPedestalWidths (0),
00018   mGains (0), mGainWidths (0),  
00019   mQieShapeCache (0), mQIEData(0), 
00020   mElectronicsMap(0),
00021   mRespCorrs(0),
00022   mL1TriggerObjects(0),
00023   mTimeCorrs(0),
00024   mLUTCorrs(0),
00025   mPFCorrs(0),
00026   mLutMetadata(0),
00027   mUpdateCalibrations (true), mUpdateCalibWidths(true)
00028  {}
00029 
00030 const HcalTopology* HcalDbService::getTopologyUsed() const {
00031   if (mPedestals && mPedestals->topo()) return mPedestals->topo();
00032   if (mGains && mGains->topo()) return mGains->topo();
00033   if (mRespCorrs && mRespCorrs->topo()) return mRespCorrs->topo();
00034   if (mL1TriggerObjects && mL1TriggerObjects->topo()) return mL1TriggerObjects->topo();
00035   if (mLutMetadata && mLutMetadata->topo()) return mLutMetadata->topo();
00036   return 0;
00037 }
00038 
00039 
00040 const HcalCalibrations& HcalDbService::getHcalCalibrations(const HcalGenericDetId& fId) const 
00041 { 
00042   if (mUpdateCalibrations) 
00043     {
00044       buildCalibrations();
00045       mUpdateCalibrations = false;
00046     }
00047   return mCalibSet.getCalibrations(fId); 
00048 }
00049 
00050 const HcalCalibrationWidths& HcalDbService::getHcalCalibrationWidths(const HcalGenericDetId& fId) const 
00051 { 
00052   if (mUpdateCalibWidths) 
00053     {
00054       buildCalibWidths();
00055       mUpdateCalibWidths = false;
00056     }
00057   return mCalibWidthSet.getCalibrationWidths(fId); 
00058 }
00059 
00060 void HcalDbService::buildCalibrations() const {
00061   // we use the set of ids for pedestals as the master list
00062   if ((!mPedestals) || (!mGains) || (!mQIEData) || (!mRespCorrs) || (!mTimeCorrs) || (!mLUTCorrs) ) return;
00063 
00064   std::vector<DetId> ids=mPedestals->getAllChannels();
00065   bool pedsInADC = mPedestals->isADC();
00066   // clear the calibrations set
00067   mCalibSet.clear();
00068   // loop!
00069   HcalCalibrations tool;
00070 
00071   //  std::cout << " length of id-vector: " << ids.size() << std::endl;
00072   for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
00073     // make
00074     bool ok=makeHcalCalibration(*id,&tool,pedsInADC);
00075     // store
00076     if (ok) mCalibSet.setCalibrations(*id,tool);
00077     //    std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
00078   }
00079   mCalibSet.sort();
00080 }
00081 
00082 void HcalDbService::buildCalibWidths() const {
00083   // we use the set of ids for pedestal widths as the master list
00084   if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData) ) return;
00085 
00086   std::vector<DetId> ids=mPedestalWidths->getAllChannels();
00087   bool pedsInADC = mPedestalWidths->isADC();
00088   // clear the calibrations set
00089   mCalibWidthSet.clear();
00090   // loop!
00091   HcalCalibrationWidths tool;
00092 
00093   //  std::cout << " length of id-vector: " << ids.size() << std::endl;
00094   for (std::vector<DetId>::const_iterator id=ids.begin(); id!=ids.end(); ++id) {
00095     // make
00096     bool ok=makeHcalCalibrationWidth(*id,&tool,pedsInADC);
00097     // store
00098     if (ok) mCalibWidthSet.setCalibrationWidths(*id,tool);
00099     //    std::cout << "Hcal calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
00100   }
00101   mCalibWidthSet.sort();
00102 }
00103 
00104 bool HcalDbService::makeHcalCalibration (const HcalGenericDetId& fId, HcalCalibrations* fObject, bool pedestalInADC) const {
00105   if (fObject) {
00106     const HcalPedestal* pedestal = getPedestal (fId);
00107     const HcalGain* gain = getGain (fId);
00108     const HcalRespCorr* respcorr = getHcalRespCorr (fId);
00109     const HcalTimeCorr* timecorr = getHcalTimeCorr (fId);
00110     const HcalLUTCorr* lutcorr = getHcalLUTCorr (fId);
00111 
00112     if (pedestalInADC) {
00113       const HcalQIECoder* coder=getHcalCoder(fId);
00114       const HcalQIEShape* shape=getHcalShape(coder);
00115       if (pedestal && gain && shape && coder && respcorr && timecorr && lutcorr) {
00116         float pedTrue[4];
00117         for (int i=0; i<4; i++) {
00118           float x=pedestal->getValues()[i];
00119           int x1=(int)std::floor(x);
00120           int x2=(int)std::floor(x+1);
00121           // y = (y2-y1)/(x2-x1) * (x - x1) + y1  [note: x2-x1=1]
00122           float y2=coder->charge(*shape,x2,i);
00123           float y1=coder->charge(*shape,x1,i);
00124           pedTrue[i]=(y2-y1)*(x-x1)+y1;
00125         }
00126         *fObject = HcalCalibrations (gain->getValues (), pedTrue, respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
00127         return true; 
00128       }
00129     } else {
00130       if (pedestal && gain && respcorr && timecorr && lutcorr) {
00131         *fObject = HcalCalibrations (gain->getValues (), pedestal->getValues (), respcorr->getValue(), timecorr->getValue(), lutcorr->getValue() );
00132         return true;
00133       }
00134     }
00135   }
00136   return false;
00137 }
00138 
00139 bool HcalDbService::makeHcalCalibrationWidth (const HcalGenericDetId& fId, 
00140                                               HcalCalibrationWidths* fObject, bool pedestalInADC) const {
00141   if (fObject) {
00142     const HcalPedestalWidth* pedestalwidth = getPedestalWidth (fId);
00143     const HcalGainWidth* gainwidth = getGainWidth (fId);
00144     if (pedestalInADC) {
00145       const HcalQIECoder* coder=getHcalCoder(fId);
00146       const HcalQIEShape* shape=getHcalShape(coder);
00147       if (pedestalwidth && gainwidth && shape && coder) {
00148         float pedTrueWidth[4];
00149         for (int i=0; i<4; i++) {
00150           float x=pedestalwidth->getWidth(i);
00151           // assume QIE is linear in low range and use x1=0 and x2=1
00152           // y = (y2-y1) * (x) [do not add any constant, only scale!]
00153           float y2=coder->charge(*shape,1,i);
00154           float y1=coder->charge(*shape,0,i);
00155           pedTrueWidth[i]=(y2-y1)*x;
00156         }
00157         *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedTrueWidth);
00158         return true; 
00159       } 
00160     } else {
00161       if (pedestalwidth && gainwidth) {
00162         float pedestalWidth [4];
00163         for (int i = 0; i < 4; i++) pedestalWidth [i] = pedestalwidth->getWidth (i);
00164         *fObject = HcalCalibrationWidths (gainwidth->getValues (), pedestalWidth);
00165         return true;
00166       }      
00167     }
00168   }
00169   return false;
00170 }  
00171 
00172 
00173 const HcalRespCorr* HcalDbService::getHcalRespCorr (const HcalGenericDetId& fId) const {
00174   if (mRespCorrs) {
00175     return mRespCorrs->getValues (fId);
00176   }
00177   return 0;
00178 }
00179 
00180 const HcalPedestal* HcalDbService::getPedestal (const HcalGenericDetId& fId) const {
00181   if (mPedestals) {
00182     return mPedestals->getValues (fId);
00183   }
00184   return 0;
00185 }
00186 
00187   const HcalPedestalWidth* HcalDbService::getPedestalWidth (const HcalGenericDetId& fId) const {
00188   if (mPedestalWidths) {
00189     return mPedestalWidths->getValues (fId);
00190   }
00191   return 0;
00192 }
00193 
00194 const HcalGain* HcalDbService::getGain (const HcalGenericDetId& fId) const {
00195   if (mGains) {
00196     return mGains->getValues(fId);
00197   }
00198   return 0;
00199 }
00200 
00201   const HcalGainWidth* HcalDbService::getGainWidth (const HcalGenericDetId& fId) const {
00202   if (mGainWidths) {
00203     return mGainWidths->getValues (fId);
00204   }
00205   return 0;
00206 }
00207 
00208 const HcalQIECoder* HcalDbService::getHcalCoder (const HcalGenericDetId& fId) const {
00209   if (mQIEData) {
00210     return mQIEData->getCoder (fId);
00211   }
00212   return 0;
00213 }
00214 
00215 const HcalQIEShape* HcalDbService::getHcalShape (const HcalGenericDetId& fId) const {
00216   if (mQIEData) {
00217     return &mQIEData->getShape (fId);
00218   }
00219   return 0;
00220 }
00221 
00222 const HcalQIEShape* HcalDbService::getHcalShape (const HcalQIECoder *coder) const {
00223   if (mQIEData) {
00224     return &mQIEData->getShape(coder);
00225   }
00226   return 0;
00227 }
00228 
00229 
00230 const HcalElectronicsMap* HcalDbService::getHcalMapping () const {
00231   return mElectronicsMap;
00232 }
00233 
00234 const HcalL1TriggerObject* HcalDbService::getHcalL1TriggerObject (const HcalGenericDetId& fId) const
00235 {
00236   return mL1TriggerObjects->getValues (fId);
00237 }
00238 
00239 const HcalChannelStatus* HcalDbService::getHcalChannelStatus (const HcalGenericDetId& fId) const
00240 {
00241   return mChannelQuality->getValues (fId);
00242 }
00243 
00244 const HcalZSThreshold* HcalDbService::getHcalZSThreshold (const HcalGenericDetId& fId) const
00245 {
00246   return mZSThresholds->getValues (fId);
00247 }
00248 
00249 const HcalTimeCorr* HcalDbService::getHcalTimeCorr (const HcalGenericDetId& fId) const {
00250   if (mTimeCorrs) {
00251     return mTimeCorrs->getValues (fId);
00252   }
00253   return 0;
00254 }
00255 
00256 const HcalLUTCorr* HcalDbService::getHcalLUTCorr (const HcalGenericDetId& fId) const {
00257   if (mLUTCorrs) {
00258     return mLUTCorrs->getValues (fId);
00259   }
00260   return 0;
00261 }
00262 
00263 const HcalPFCorr* HcalDbService::getHcalPFCorr (const HcalGenericDetId& fId) const {
00264   if (mPFCorrs) {
00265     return mPFCorrs->getValues (fId);
00266   }
00267   return 0;
00268 }
00269 
00270 const HcalLutMetadata* HcalDbService::getHcalLutMetadata () const {
00271   return mLutMetadata;
00272 }
00273 
00274 TYPELOOKUP_DATA_REG(HcalDbService);