Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 //
0003 // F.Ratnikov (UMd), Aug. 9, 2005
0004 // Adapted for CASTOR by L. Mundim
0005 //
0006 
0007 #ifndef CastorDbService_h
0008 #define CastorDbService_h
0009 
0010 #include <memory>
0011 #include <map>
0012 
0013 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0014 #include "CalibFormats/CastorObjects/interface/CastorChannelCoder.h"
0015 #include "CalibFormats/CastorObjects/interface/QieShape.h"
0016 #include "CalibFormats/CastorObjects/interface/CastorCoder.h"
0017 #include "CalibFormats/CastorObjects/interface/CastorCalibrationsSet.h"
0018 #include "CalibFormats/CastorObjects/interface/CastorCalibrationWidthsSet.h"
0019 
0020 #include "CondFormats/CastorObjects/interface/AllObjects.h"
0021 
0022 class CastorCalibrations;
0023 class CastorCalibrationWidths;
0024 
0025 class CastorDbService {
0026 public:
0027   CastorDbService();
0028 
0029   const CastorCalibrations& getCastorCalibrations(const HcalGenericDetId& fId) const {
0030     return mCalibSet.getCalibrations(fId);
0031   }
0032   const CastorCalibrationWidths& getCastorCalibrationWidths(const HcalGenericDetId& fId) const {
0033     return mCalibWidthSet.getCalibrationWidths(fId);
0034   }
0035 
0036   const CastorPedestal* getPedestal(const HcalGenericDetId& fId) const;
0037   const CastorPedestalWidth* getPedestalWidth(const HcalGenericDetId& fId) const;
0038   const CastorGain* getGain(const HcalGenericDetId& fId) const;
0039   const CastorGainWidth* getGainWidth(const HcalGenericDetId& fId) const;
0040   const CastorQIECoder* getCastorCoder(const HcalGenericDetId& fId) const;
0041   const CastorQIEShape* getCastorShape() const;
0042   const CastorElectronicsMap* getCastorMapping() const;
0043   const CastorChannelStatus* getCastorChannelStatus(const HcalGenericDetId& fId) const;
0044 
0045   void setData(const CastorPedestals* fItem) { mPedestals = fItem; }
0046   void setData(const CastorPedestalWidths* fItem) { mPedestalWidths = fItem; }
0047   void setData(const CastorGains* fItem) { mGains = fItem; }
0048   void setData(const CastorGainWidths* fItem) { mGainWidths = fItem; }
0049   void setData(const CastorQIEData* fItem) { mQIEData = fItem; }
0050   void setData(const CastorChannelQuality* fItem) { mChannelQuality = fItem; }
0051   void setData(const CastorElectronicsMap* fItem) { mElectronicsMap = fItem; }
0052 
0053   void buildCalibrations();
0054   void buildCalibWidths();
0055 
0056 private:
0057   bool makeCastorCalibration(const HcalGenericDetId& fId, CastorCalibrations* fObject, bool pedestalInADC) const;
0058   bool makeCastorCalibrationWidth(const HcalGenericDetId& fId,
0059                                   CastorCalibrationWidths* fObject,
0060                                   bool pedestalInADC) const;
0061   const CastorPedestals* mPedestals;
0062   const CastorPedestalWidths* mPedestalWidths;
0063   const CastorGains* mGains;
0064   const CastorGainWidths* mGainWidths;
0065   const CastorQIEData* mQIEData;
0066   const CastorChannelQuality* mChannelQuality;
0067   const CastorElectronicsMap* mElectronicsMap;
0068   CastorCalibrationsSet mCalibSet;
0069   CastorCalibrationWidthsSet mCalibWidthSet;
0070 };
0071 
0072 #endif