Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include <memory>
00008 #include <iostream>
00009
00010 #include "FWCore/Framework/interface/ValidityInterval.h"
00011 #include "DataFormats/HcalDetId/interface/HcalCastorDetId.h"
00012 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00013 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
00014 #include "CalibCalorimetry/CastorCalib/interface/CastorDbHardcode.h"
00015
00016 #include "CondFormats/DataRecord/interface/CastorPedestalsRcd.h"
00017 #include "CondFormats/DataRecord/interface/CastorPedestalWidthsRcd.h"
00018 #include "CondFormats/DataRecord/interface/CastorGainsRcd.h"
00019 #include "CondFormats/DataRecord/interface/CastorGainWidthsRcd.h"
00020 #include "CondFormats/DataRecord/interface/CastorElectronicsMapRcd.h"
00021 #include "CondFormats/DataRecord/interface/CastorChannelQualityRcd.h"
00022 #include "CondFormats/DataRecord/interface/CastorQIEDataRcd.h"
00023 #include "CondFormats/DataRecord/interface/CastorRecoParamsRcd.h"
00024 #include "CondFormats/DataRecord/interface/CastorSaturationCorrsRcd.h"
00025
00026 #include "Geometry/ForwardGeometry/interface/CastorTopology.h"
00027 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029
00030 #include "CastorHardcodeCalibrations.h"
00031
00032
00033
00034
00035 using namespace cms;
00036
00037 namespace {
00038
00039 std::vector<HcalGenericDetId> allCells (bool h2_mode) {
00040 static std::vector<HcalGenericDetId> result;
00041 if (result.size () <= 0) {
00042
00043 CastorTopology castortopology;
00044 HcalCastorDetId cell;
00045 HcalCastorDetId::Section section = HcalCastorDetId::EM;
00046
00047 for(int sector=1; sector<17; sector++) {
00048 for(int module=1; module<3; module++) {
00049 cell = HcalCastorDetId(section, true, sector, module);
00050 if (castortopology.valid(cell)) result.push_back(cell);
00051 cell = HcalCastorDetId(section, false, sector, module);
00052 if (castortopology.valid(cell)) result.push_back(cell);
00053 }
00054 }
00055
00056 section = HcalCastorDetId::HAD;
00057 for(int sector= 1; sector < 17; sector++){
00058 for(int module=3; module<15; module++) {
00059 cell = HcalCastorDetId(section, true, sector, module);
00060 if(castortopology.valid(cell)) result.push_back(cell);
00061 cell = HcalCastorDetId(section, false, sector, module);
00062 if(castortopology.valid(cell)) result.push_back(cell);
00063 }
00064 }
00065
00066 }
00067 return result;
00068
00069 }
00070 }
00071
00072 CastorHardcodeCalibrations::CastorHardcodeCalibrations ( const edm::ParameterSet& iConfig )
00073
00074 {
00075 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::CastorHardcodeCalibrations->...";
00076
00077 h2mode_=iConfig.getUntrackedParameter<bool>("H2Mode",false);
00078 std::vector <std::string> toGet = iConfig.getUntrackedParameter <std::vector <std::string> > ("toGet");
00079 for(std::vector <std::string>::iterator objectName = toGet.begin(); objectName != toGet.end(); ++objectName ) {
00080 bool all = *objectName == "all";
00081 if ((*objectName == "Pedestals") || all) {
00082 setWhatProduced (this, &CastorHardcodeCalibrations::producePedestals);
00083 findingRecord <CastorPedestalsRcd> ();
00084 }
00085 if ((*objectName == "PedestalWidths") || all) {
00086 setWhatProduced (this, &CastorHardcodeCalibrations::producePedestalWidths);
00087 findingRecord <CastorPedestalWidthsRcd> ();
00088 }
00089 if ((*objectName == "Gains") || all) {
00090 setWhatProduced (this, &CastorHardcodeCalibrations::produceGains);
00091 findingRecord <CastorGainsRcd> ();
00092 }
00093 if ((*objectName == "GainWidths") || all) {
00094 setWhatProduced (this, &CastorHardcodeCalibrations::produceGainWidths);
00095 findingRecord <CastorGainWidthsRcd> ();
00096 }
00097 if ((*objectName == "QIEData") || all) {
00098 setWhatProduced (this, &CastorHardcodeCalibrations::produceQIEData);
00099 findingRecord <CastorQIEDataRcd> ();
00100 }
00101 if ((*objectName == "ChannelQuality") || (*objectName == "channelQuality") || all) {
00102 setWhatProduced (this, &CastorHardcodeCalibrations::produceChannelQuality);
00103 findingRecord <CastorChannelQualityRcd> ();
00104 }
00105 if ((*objectName == "ElectronicsMap") || (*objectName == "electronicsMap") || all) {
00106 setWhatProduced (this, &CastorHardcodeCalibrations::produceElectronicsMap);
00107 findingRecord <CastorElectronicsMapRcd> ();
00108 }
00109 if ((*objectName == "RecoParams") || all) {
00110 setWhatProduced (this, &CastorHardcodeCalibrations::produceRecoParams);
00111 findingRecord <CastorRecoParamsRcd> ();
00112 }
00113 if ((*objectName == "SaturationCorrs") || all) {
00114 setWhatProduced (this, &CastorHardcodeCalibrations::produceSaturationCorrs);
00115 findingRecord <CastorSaturationCorrsRcd> ();
00116 }
00117 }
00118 }
00119
00120
00121 CastorHardcodeCalibrations::~CastorHardcodeCalibrations()
00122 {
00123 }
00124
00125
00126
00127
00128
00129 void
00130 CastorHardcodeCalibrations::setIntervalFor( const edm::eventsetup::EventSetupRecordKey& iKey, const edm::IOVSyncValue& iTime, edm::ValidityInterval& oInterval ) {
00131 std::string record = iKey.name ();
00132 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::setIntervalFor-> key: " << record << " time: " << iTime.eventID() << '/' << iTime.time ().value ();
00133 oInterval = edm::ValidityInterval (edm::IOVSyncValue::beginOfTime(), edm::IOVSyncValue::endOfTime());
00134 }
00135
00136 std::auto_ptr<CastorPedestals> CastorHardcodeCalibrations::producePedestals (const CastorPedestalsRcd&) {
00137 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::producePedestals-> ...";
00138 std::auto_ptr<CastorPedestals> result (new CastorPedestals (false));
00139 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00140 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00141 CastorPedestal item = CastorDbHardcode::makePedestal (*cell);
00142 result->addValues(item);
00143 }
00144 return result;
00145 }
00146
00147 std::auto_ptr<CastorPedestalWidths> CastorHardcodeCalibrations::producePedestalWidths (const CastorPedestalWidthsRcd&) {
00148 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::producePedestalWidths-> ...";
00149 std::auto_ptr<CastorPedestalWidths> result (new CastorPedestalWidths (false));
00150 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00151 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00152 CastorPedestalWidth item = CastorDbHardcode::makePedestalWidth (*cell);
00153 result->addValues(item);
00154 }
00155 return result;
00156 }
00157
00158 std::auto_ptr<CastorGains> CastorHardcodeCalibrations::produceGains (const CastorGainsRcd&) {
00159 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceGains-> ...";
00160 std::auto_ptr<CastorGains> result (new CastorGains ());
00161 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00162 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00163 CastorGain item = CastorDbHardcode::makeGain (*cell);
00164 result->addValues(item);
00165 }
00166 return result;
00167 }
00168
00169 std::auto_ptr<CastorGainWidths> CastorHardcodeCalibrations::produceGainWidths (const CastorGainWidthsRcd&) {
00170 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceGainWidths-> ...";
00171 std::auto_ptr<CastorGainWidths> result (new CastorGainWidths ());
00172 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00173 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00174 CastorGainWidth item = CastorDbHardcode::makeGainWidth (*cell);
00175 result->addValues(item);
00176 }
00177 return result;
00178 }
00179
00180 std::auto_ptr<CastorQIEData> CastorHardcodeCalibrations::produceQIEData (const CastorQIEDataRcd& rcd) {
00181 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceQIEData-> ...";
00182 std::auto_ptr<CastorQIEData> result (new CastorQIEData ());
00183 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00184 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00185 CastorQIECoder coder = CastorDbHardcode::makeQIECoder (*cell);
00186 result->addCoder(coder);
00187 }
00188 return result;
00189 }
00190
00191 std::auto_ptr<CastorChannelQuality> CastorHardcodeCalibrations::produceChannelQuality (const CastorChannelQualityRcd& rcd) {
00192 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceChannelQuality-> ...";
00193 std::auto_ptr<CastorChannelQuality> result (new CastorChannelQuality ());
00194 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00195 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00196 CastorChannelStatus item(cell->rawId(),CastorChannelStatus::GOOD);
00197 result->addValues(item);
00198 }
00199 return result;
00200 }
00201
00202 std::auto_ptr<CastorElectronicsMap> CastorHardcodeCalibrations::produceElectronicsMap (const CastorElectronicsMapRcd& rcd) {
00203 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceElectronicsMap-> ...";
00204
00205 std::auto_ptr<CastorElectronicsMap> result (new CastorElectronicsMap ());
00206 CastorDbHardcode::makeHardcodeMap(*result);
00207 return result;
00208 }
00209
00210 std::auto_ptr<CastorRecoParams> CastorHardcodeCalibrations::produceRecoParams (const CastorRecoParamsRcd& rcd) {
00211 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceRecoParams-> ...";
00212 std::auto_ptr<CastorRecoParams> result (new CastorRecoParams ());
00213 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00214 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00215 CastorRecoParam item = CastorDbHardcode::makeRecoParam (*cell);
00216 result->addValues(item);
00217 }
00218 return result;
00219 }
00220
00221 std::auto_ptr<CastorSaturationCorrs> CastorHardcodeCalibrations::produceSaturationCorrs (const CastorSaturationCorrsRcd& rcd) {
00222 edm::LogInfo("HCAL") << "CastorHardcodeCalibrations::produceSaturationCorrs-> ...";
00223 std::auto_ptr<CastorSaturationCorrs> result (new CastorSaturationCorrs ());
00224 std::vector <HcalGenericDetId> cells = allCells(h2mode_);
00225 for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00226 CastorSaturationCorr item = CastorDbHardcode::makeSaturationCorr (*cell);
00227 result->addValues(item);
00228 }
00229 return result;
00230 }