Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:11

0001 //
0002 // F.Ratnikov (UMd), Aug. 9, 2005
0003 // Adapted for CASTOR by L. Mundim
0004 //
0005 
0006 #include "FWCore/Utilities/interface/typelookup.h"
0007 
0008 #include "CalibFormats/CastorObjects/interface/CastorDbService.h"
0009 #include "CalibFormats/CastorObjects/interface/CastorCoderDb.h"
0010 #include "CalibFormats/CastorObjects/interface/QieShape.h"
0011 
0012 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0013 
0014 #include <cmath>
0015 
0016 CastorDbService::CastorDbService()
0017     : mPedestals(nullptr),
0018       mPedestalWidths(nullptr),
0019       mGains(nullptr),
0020       mGainWidths(nullptr),
0021       mQIEData(nullptr),
0022       mChannelQuality(nullptr),
0023       mElectronicsMap(nullptr) {}
0024 
0025 bool CastorDbService::makeCastorCalibration(const HcalGenericDetId& fId,
0026                                             CastorCalibrations* fObject,
0027                                             bool pedestalInADC) const {
0028   if (fObject) {
0029     const CastorPedestal* pedestal = getPedestal(fId);
0030     const CastorGain* gain = getGain(fId);
0031 
0032     if (pedestalInADC) {
0033       const CastorQIEShape* shape = getCastorShape();
0034       const CastorQIECoder* coder = getCastorCoder(fId);
0035       if (pedestal && gain && shape && coder) {
0036         float pedTrue[4];
0037         for (int i = 0; i < 4; i++) {
0038           float x = pedestal->getValues()[i];
0039           int x1 = (int)std::floor(x);
0040           int x2 = (int)std::floor(x + 1);
0041           // y = (y2-y1)/(x2-x1) * (x - x1) + y1  [note: x2-x1=1]
0042           float y2 = coder->charge(*shape, x2, i);
0043           float y1 = coder->charge(*shape, x1, i);
0044           pedTrue[i] = (y2 - y1) * (x - x1) + y1;
0045         }
0046         *fObject = CastorCalibrations(gain->getValues(), pedTrue);
0047         return true;
0048       }
0049     } else {
0050       if (pedestal && gain) {
0051         *fObject = CastorCalibrations(gain->getValues(), pedestal->getValues());
0052         return true;
0053       }
0054     }
0055   }
0056   return false;
0057 }
0058 
0059 void CastorDbService::buildCalibrations() {
0060   // we use the set of ids for pedestals as the master list
0061   if ((!mPedestals) || (!mGains) || (!mQIEData))
0062     return;
0063   std::vector<DetId> ids = mPedestals->getAllChannels();
0064   bool pedsInADC = mPedestals->isADC();
0065   // clear the calibrations set
0066   mCalibSet.clear();
0067   // loop!
0068   CastorCalibrations tool;
0069 
0070   //  std::cout << " length of id-vector: " << ids.size() << std::endl;
0071   for (std::vector<DetId>::const_iterator id = ids.begin(); id != ids.end(); ++id) {
0072     // make
0073     bool ok = makeCastorCalibration(*id, &tool, pedsInADC);
0074     // store
0075     if (ok)
0076       mCalibSet.setCalibrations(*id, tool);
0077     //    std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
0078   }
0079   mCalibSet.sort();
0080 }
0081 
0082 void CastorDbService::buildCalibWidths() {
0083   // we use the set of ids for pedestal widths as the master list
0084   if ((!mPedestalWidths) || (!mGainWidths) || (!mQIEData))
0085     return;
0086 
0087   std::vector<DetId> ids = mPedestalWidths->getAllChannels();
0088   bool pedsInADC = mPedestalWidths->isADC();
0089   // clear the calibrations set
0090   mCalibWidthSet.clear();
0091   // loop!
0092   CastorCalibrationWidths tool;
0093 
0094   //  std::cout << " length of id-vector: " << ids.size() << std::endl;
0095   for (std::vector<DetId>::const_iterator id = ids.begin(); id != ids.end(); ++id) {
0096     // make
0097     bool ok = makeCastorCalibrationWidth(*id, &tool, pedsInADC);
0098     // store
0099     if (ok)
0100       mCalibWidthSet.setCalibrationWidths(*id, tool);
0101     //    std::cout << "Castor calibrations built... detid no. " << HcalGenericDetId(*id) << std::endl;
0102   }
0103   mCalibWidthSet.sort();
0104 }
0105 
0106 bool CastorDbService::makeCastorCalibrationWidth(const HcalGenericDetId& fId,
0107                                                  CastorCalibrationWidths* fObject,
0108                                                  bool pedestalInADC) const {
0109   if (fObject) {
0110     const CastorPedestalWidth* pedestalwidth = getPedestalWidth(fId);
0111     const CastorGainWidth* gainwidth = getGainWidth(fId);
0112     if (pedestalInADC) {
0113       const CastorQIEShape* shape = getCastorShape();
0114       const CastorQIECoder* coder = getCastorCoder(fId);
0115       if (pedestalwidth && gainwidth && shape && coder) {
0116         float pedTrueWidth[4];
0117         for (int i = 0; i < 4; i++) {
0118           float x = pedestalwidth->getWidth(i);
0119           // assume QIE is linear in low range and use x1=0 and x2=1
0120           // y = (y2-y1) * (x) [do not add any constant, only scale!]
0121           float y2 = coder->charge(*shape, 1, i);
0122           float y1 = coder->charge(*shape, 0, i);
0123           pedTrueWidth[i] = (y2 - y1) * x;
0124         }
0125         *fObject = CastorCalibrationWidths(gainwidth->getValues(), pedTrueWidth);
0126         return true;
0127       }
0128     } else {
0129       if (pedestalwidth && gainwidth) {
0130         float pedestalWidth[4];
0131         for (int i = 0; i < 4; i++)
0132           pedestalWidth[i] = pedestalwidth->getWidth(i);
0133         *fObject = CastorCalibrationWidths(gainwidth->getValues(), pedestalWidth);
0134         return true;
0135       }
0136     }
0137   }
0138   return false;
0139 }
0140 
0141 const CastorPedestal* CastorDbService::getPedestal(const HcalGenericDetId& fId) const {
0142   if (mPedestals) {
0143     return mPedestals->getValues(fId);
0144   }
0145   return nullptr;
0146 }
0147 
0148 const CastorPedestalWidth* CastorDbService::getPedestalWidth(const HcalGenericDetId& fId) const {
0149   if (mPedestalWidths) {
0150     return mPedestalWidths->getValues(fId);
0151   }
0152   return nullptr;
0153 }
0154 
0155 const CastorGain* CastorDbService::getGain(const HcalGenericDetId& fId) const {
0156   if (mGains) {
0157     return mGains->getValues(fId);
0158   }
0159   return nullptr;
0160 }
0161 
0162 const CastorGainWidth* CastorDbService::getGainWidth(const HcalGenericDetId& fId) const {
0163   if (mGainWidths) {
0164     return mGainWidths->getValues(fId);
0165   }
0166   return nullptr;
0167 }
0168 
0169 const CastorQIECoder* CastorDbService::getCastorCoder(const HcalGenericDetId& fId) const {
0170   if (mQIEData) {
0171     return mQIEData->getCoder(fId);
0172   }
0173   return nullptr;
0174 }
0175 
0176 const CastorQIEShape* CastorDbService::getCastorShape() const {
0177   if (mQIEData) {
0178     return &mQIEData->getShape();
0179   }
0180   return nullptr;
0181 }
0182 
0183 const CastorElectronicsMap* CastorDbService::getCastorMapping() const { return mElectronicsMap; }
0184 
0185 const CastorChannelStatus* CastorDbService::getCastorChannelStatus(const HcalGenericDetId& fId) const {
0186   return mChannelQuality->getValues(fId);
0187 }
0188 
0189 TYPELOOKUP_DATA_REG(CastorDbService);