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