CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/CalibCalorimetry/HcalPlugins/src/HcalHardcodeCalibrations.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 // Original Author:  Fedor Ratnikov
00003 // $Id: HcalHardcodeCalibrations.cc,v 1.34 2012/11/12 20:42:39 dlange Exp $
00004 //
00005 //
00006 
00007 #include <memory>
00008 #include <iostream>
00009 
00010 #include "FWCore/Framework/interface/ValidityInterval.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 #include "DataFormats/HcalDetId/interface/HcalZDCDetId.h"
00013 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00014 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
00015 #include "CalibCalorimetry/HcalAlgos/interface/HcalDbHardcode.h"
00016 
00017 #include "CondFormats/DataRecord/interface/HcalAllRcds.h"
00018 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00019 
00020 #include "Geometry/ForwardGeometry/interface/ZdcTopology.h"
00021 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 
00024 #include "HcalHardcodeCalibrations.h"
00025 
00026 // class decleration
00027 //
00028 
00029 using namespace cms;
00030 
00031 namespace {
00032 
00033   std::vector<HcalGenericDetId> allCells (const HcalTopology& hcaltopology) {
00034   static std::vector<HcalGenericDetId> result;
00035   int maxDepthHB=hcaltopology.maxDepthHB();
00036   int maxDepthHE=hcaltopology.maxDepthHE();
00037   if (result.size () <= 0) {
00038     for (int eta = -50; eta < 50; eta++) {
00039       for (int phi = 0; phi < 100; phi++) {
00040         for (int depth = 1; depth < maxDepthHB + maxDepthHE; depth++) {
00041           for (int det = 1; det < 5; det++) {
00042             HcalDetId cell ((HcalSubdetector) det, eta, phi, depth);
00043             if (hcaltopology.valid(cell)) result.push_back (cell);
00044           }
00045         }
00046       }
00047     } 
00048     ZdcTopology zdctopology;
00049     HcalZDCDetId zcell;
00050     HcalZDCDetId::Section section  = HcalZDCDetId::EM;
00051     for(int depth= 1; depth < 6; depth++){
00052       zcell = HcalZDCDetId(section, true, depth);
00053       if(zdctopology.valid(zcell)) result.push_back(zcell);
00054       zcell = HcalZDCDetId(section, false, depth);
00055       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00056     }
00057     section = HcalZDCDetId::HAD;
00058     for(int depth= 1; depth < 5; depth++){
00059       zcell = HcalZDCDetId(section, true, depth);
00060       if(zdctopology.valid(zcell)) result.push_back(zcell);
00061       zcell = HcalZDCDetId(section, false, depth);
00062       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00063     }
00064     section = HcalZDCDetId::LUM;
00065     for(int depth= 1; depth < 3; depth++){
00066       zcell = HcalZDCDetId(section, true, depth);
00067       if(zdctopology.valid(zcell)) result.push_back(zcell);
00068       zcell = HcalZDCDetId(section, false, depth);
00069       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00070     }
00071   }
00072   return result;
00073 }
00074 
00075 }
00076 
00077 HcalHardcodeCalibrations::HcalHardcodeCalibrations ( const edm::ParameterSet& iConfig )
00078 {
00079   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::HcalHardcodeCalibrations->...";
00080   
00081   std::vector <std::string> toGet = iConfig.getUntrackedParameter <std::vector <std::string> > ("toGet");
00082   for(std::vector <std::string>::iterator objectName = toGet.begin(); objectName != toGet.end(); ++objectName ) {
00083     bool all = *objectName == "all";
00084     if ((*objectName == "Pedestals") || all) {
00085       setWhatProduced (this, &HcalHardcodeCalibrations::producePedestals);
00086       findingRecord <HcalPedestalsRcd> ();
00087     }
00088     if ((*objectName == "PedestalWidths") || all) {
00089       setWhatProduced (this, &HcalHardcodeCalibrations::producePedestalWidths);
00090       findingRecord <HcalPedestalWidthsRcd> ();
00091     }
00092     if ((*objectName == "Gains") || all) {
00093       setWhatProduced (this, &HcalHardcodeCalibrations::produceGains);
00094       findingRecord <HcalGainsRcd> ();
00095     }
00096     if ((*objectName == "GainWidths") || all) {
00097       setWhatProduced (this, &HcalHardcodeCalibrations::produceGainWidths);
00098       findingRecord <HcalGainWidthsRcd> ();
00099     }
00100     if ((*objectName == "QIEData") || all) {
00101       setWhatProduced (this, &HcalHardcodeCalibrations::produceQIEData);
00102       findingRecord <HcalQIEDataRcd> ();
00103     }
00104     if ((*objectName == "ChannelQuality") || (*objectName == "channelQuality") || all) {
00105       setWhatProduced (this, &HcalHardcodeCalibrations::produceChannelQuality);
00106       findingRecord <HcalChannelQualityRcd> ();
00107     }
00108     if ((*objectName == "ElectronicsMap") || (*objectName == "electronicsMap") || all) {
00109       setWhatProduced (this, &HcalHardcodeCalibrations::produceElectronicsMap);
00110       findingRecord <HcalElectronicsMapRcd> ();
00111     }
00112     if ((*objectName == "ZSThresholds") || (*objectName == "zsThresholds") || all) {
00113       setWhatProduced (this, &HcalHardcodeCalibrations::produceZSThresholds);
00114       findingRecord <HcalZSThresholdsRcd> ();
00115     }
00116     if ((*objectName == "RespCorrs") || (*objectName == "ResponseCorrection") || all) {
00117       setWhatProduced (this, &HcalHardcodeCalibrations::produceRespCorrs);
00118       findingRecord <HcalRespCorrsRcd> ();
00119     }
00120     if ((*objectName == "LUTCorrs") || (*objectName == "LUTCorrection") || all) {
00121       setWhatProduced (this, &HcalHardcodeCalibrations::produceLUTCorrs);
00122       findingRecord <HcalLUTCorrsRcd> ();
00123     }
00124     if ((*objectName == "PFCorrs") || (*objectName == "PFCorrection") || all) {
00125       setWhatProduced (this, &HcalHardcodeCalibrations::producePFCorrs);
00126       findingRecord <HcalPFCorrsRcd> ();
00127     }
00128     if ((*objectName == "TimeCorrs") || (*objectName == "TimeCorrection") || all) {
00129       setWhatProduced (this, &HcalHardcodeCalibrations::produceTimeCorrs);
00130       findingRecord <HcalTimeCorrsRcd> ();
00131     }
00132     if ((*objectName == "L1TriggerObjects") || (*objectName == "L1Trigger") || all) {
00133       setWhatProduced (this, &HcalHardcodeCalibrations::produceL1TriggerObjects);
00134       findingRecord <HcalL1TriggerObjectsRcd> ();
00135     }
00136     if ((*objectName == "ValidationCorrs") || (*objectName == "ValidationCorrection") || all) {
00137       setWhatProduced (this, &HcalHardcodeCalibrations::produceValidationCorrs);
00138       findingRecord <HcalValidationCorrsRcd> ();
00139     }
00140     if ((*objectName == "LutMetadata") || (*objectName == "lutMetadata") || all) {
00141       setWhatProduced (this, &HcalHardcodeCalibrations::produceLutMetadata);
00142       findingRecord <HcalLutMetadataRcd> ();
00143     }
00144     if ((*objectName == "DcsValues") || all) {
00145       setWhatProduced (this, &HcalHardcodeCalibrations::produceDcsValues);
00146       findingRecord <HcalDcsRcd> ();
00147     }
00148     if ((*objectName == "DcsMap") || (*objectName == "dcsMap") || all) {
00149       setWhatProduced (this, &HcalHardcodeCalibrations::produceDcsMap);
00150       findingRecord <HcalDcsMapRcd> ();
00151     }
00152     if ((*objectName == "RecoParams") || all) {
00153       setWhatProduced (this, &HcalHardcodeCalibrations::produceRecoParams);
00154       findingRecord <HcalRecoParamsRcd> ();
00155     }
00156     if ((*objectName == "LongRecoParams") || all) {
00157       setWhatProduced (this, &HcalHardcodeCalibrations::produceLongRecoParams);
00158       findingRecord <HcalLongRecoParamsRcd> ();
00159     }
00160     if ((*objectName == "MCParams") || all) {
00161       setWhatProduced (this, &HcalHardcodeCalibrations::produceMCParams);
00162       findingRecord <HcalMCParamsRcd> ();
00163     }
00164     if ((*objectName == "FlagHFDigiTimeParams") || all) {
00165       setWhatProduced (this, &HcalHardcodeCalibrations::produceFlagHFDigiTimeParams);
00166       findingRecord <HcalFlagHFDigiTimeParamsRcd> ();
00167     }
00168   }
00169 }
00170 
00171 
00172 HcalHardcodeCalibrations::~HcalHardcodeCalibrations()
00173 {
00174 }
00175 
00176 //
00177 // member functions
00178 //
00179 void 
00180 HcalHardcodeCalibrations::setIntervalFor( const edm::eventsetup::EventSetupRecordKey& iKey, const edm::IOVSyncValue& iTime, edm::ValidityInterval& oInterval ) {
00181   std::string record = iKey.name ();
00182   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::setIntervalFor-> key: " << record << " time: " << iTime.eventID() << '/' << iTime.time ().value ();
00183   oInterval = edm::ValidityInterval (edm::IOVSyncValue::beginOfTime(), edm::IOVSyncValue::endOfTime()); //infinite
00184 }
00185 
00186 std::auto_ptr<HcalPedestals> HcalHardcodeCalibrations::producePedestals (const HcalPedestalsRcd& rec) {
00187   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePedestals-> ...";
00188   edm::ESHandle<HcalTopology> htopo;
00189   rec.getRecord<IdealGeometryRecord>().get(htopo);
00190   const HcalTopology* topo=&(*htopo);
00191 
00192   std::auto_ptr<HcalPedestals> result (new HcalPedestals (topo,false));
00193   std::vector <HcalGenericDetId> cells = allCells(*topo);
00194   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00195     HcalPedestal item = HcalDbHardcode::makePedestal (*cell);
00196     result->addValues(item);
00197   }
00198   return result;
00199 }
00200 
00201 std::auto_ptr<HcalPedestalWidths> HcalHardcodeCalibrations::producePedestalWidths (const HcalPedestalWidthsRcd& rec) {
00202   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePedestalWidths-> ...";
00203   edm::ESHandle<HcalTopology> htopo;
00204   rec.getRecord<IdealGeometryRecord>().get(htopo);
00205   const HcalTopology* topo=&(*htopo);
00206 
00207   std::auto_ptr<HcalPedestalWidths> result (new HcalPedestalWidths (topo,false));
00208   std::vector <HcalGenericDetId> cells = allCells(*htopo);
00209   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00210     HcalPedestalWidth item = HcalDbHardcode::makePedestalWidth (*cell);
00211     result->addValues(item);
00212   }
00213   return result;
00214 }
00215 
00216 std::auto_ptr<HcalGains> HcalHardcodeCalibrations::produceGains (const HcalGainsRcd& rec) {
00217   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGains-> ...";
00218   edm::ESHandle<HcalTopology> htopo;
00219   rec.getRecord<IdealGeometryRecord>().get(htopo);
00220   const HcalTopology* topo=&(*htopo);
00221 
00222   std::auto_ptr<HcalGains> result (new HcalGains (topo));
00223   std::vector <HcalGenericDetId> cells = allCells(*topo);
00224   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00225     HcalGain item = HcalDbHardcode::makeGain (*cell);
00226     result->addValues(item);
00227   }
00228   return result;
00229 }
00230 
00231 std::auto_ptr<HcalGainWidths> HcalHardcodeCalibrations::produceGainWidths (const HcalGainWidthsRcd& rec) {
00232   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGainWidths-> ...";
00233   edm::ESHandle<HcalTopology> htopo;
00234   rec.getRecord<IdealGeometryRecord>().get(htopo);
00235   const HcalTopology* topo=&(*htopo);
00236 
00237   std::auto_ptr<HcalGainWidths> result (new HcalGainWidths (topo));
00238   std::vector <HcalGenericDetId> cells = allCells(*topo);
00239   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00240     HcalGainWidth item = HcalDbHardcode::makeGainWidth (*cell);
00241     result->addValues(item);
00242   }
00243   return result;
00244 }
00245 
00246 std::auto_ptr<HcalQIEData> HcalHardcodeCalibrations::produceQIEData (const HcalQIEDataRcd& rcd) {
00247   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceQIEData-> ...";
00248   edm::ESHandle<HcalTopology> htopo;
00249   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00250   const HcalTopology* topo=&(*htopo);
00251 
00252   std::auto_ptr<HcalQIEData> result (new HcalQIEData (topo));
00253   std::vector <HcalGenericDetId> cells = allCells(*topo);
00254   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00255     HcalQIECoder coder = HcalDbHardcode::makeQIECoder (*cell);
00256     result->addCoder (coder);
00257   }
00258   return result;
00259 }
00260 
00261 std::auto_ptr<HcalChannelQuality> HcalHardcodeCalibrations::produceChannelQuality (const HcalChannelQualityRcd& rcd) {
00262   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceChannelQuality-> ...";
00263   edm::ESHandle<HcalTopology> htopo;
00264   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00265   const HcalTopology* topo=&(*htopo);
00266 
00267   std::auto_ptr<HcalChannelQuality> result (new HcalChannelQuality (topo));
00268   std::vector <HcalGenericDetId> cells = allCells(*topo);
00269   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00270     HcalChannelStatus item(cell->rawId(),0);
00271     result->addValues(item);
00272   }
00273   return result;
00274 }
00275 
00276 
00277 std::auto_ptr<HcalRespCorrs> HcalHardcodeCalibrations::produceRespCorrs (const HcalRespCorrsRcd& rcd) {
00278   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRespCorrs-> ...";
00279   edm::ESHandle<HcalTopology> htopo;
00280   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00281   const HcalTopology* topo=&(*htopo);
00282 
00283   std::auto_ptr<HcalRespCorrs> result (new HcalRespCorrs (topo));
00284   std::vector <HcalGenericDetId> cells = allCells(*topo);
00285   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00286     HcalRespCorr item(cell->rawId(),1.0);
00287     result->addValues(item);
00288   }
00289   return result;
00290 }
00291 
00292 std::auto_ptr<HcalLUTCorrs> HcalHardcodeCalibrations::produceLUTCorrs (const HcalLUTCorrsRcd& rcd) {
00293   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLUTCorrs-> ...";
00294   edm::ESHandle<HcalTopology> htopo;
00295   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00296   const HcalTopology* topo=&(*htopo);
00297 
00298   std::auto_ptr<HcalLUTCorrs> result (new HcalLUTCorrs (topo));
00299   std::vector <HcalGenericDetId> cells = allCells(*topo);
00300   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00301     HcalLUTCorr item(cell->rawId(),1.0);
00302     result->addValues(item);
00303   }
00304   return result;
00305 }
00306 
00307 std::auto_ptr<HcalPFCorrs> HcalHardcodeCalibrations::producePFCorrs (const HcalPFCorrsRcd& rcd) {
00308   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePFCorrs-> ...";
00309   edm::ESHandle<HcalTopology> htopo;
00310   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00311   const HcalTopology* topo=&(*htopo);
00312 
00313   std::auto_ptr<HcalPFCorrs> result (new HcalPFCorrs (topo));
00314   std::vector <HcalGenericDetId> cells = allCells(*topo);
00315   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00316     HcalPFCorr item(cell->rawId(),1.0);
00317     result->addValues(item);
00318   }
00319   return result;
00320 }
00321 
00322 std::auto_ptr<HcalTimeCorrs> HcalHardcodeCalibrations::produceTimeCorrs (const HcalTimeCorrsRcd& rcd) {
00323   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimeCorrs-> ...";
00324   edm::ESHandle<HcalTopology> htopo;
00325   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00326   const HcalTopology* topo=&(*htopo);
00327 
00328   std::auto_ptr<HcalTimeCorrs> result (new HcalTimeCorrs (topo));
00329   std::vector <HcalGenericDetId> cells = allCells(*topo);
00330   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00331     HcalTimeCorr item(cell->rawId(),0.0);
00332     result->addValues(item);
00333   }
00334   return result;
00335 }
00336 
00337 std::auto_ptr<HcalZSThresholds> HcalHardcodeCalibrations::produceZSThresholds (const HcalZSThresholdsRcd& rcd) {
00338   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceZSThresholds-> ...";
00339   edm::ESHandle<HcalTopology> htopo;
00340   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00341   const HcalTopology* topo=&(*htopo);
00342 
00343   std::auto_ptr<HcalZSThresholds> result (new HcalZSThresholds (topo));
00344   std::vector <HcalGenericDetId> cells = allCells(*topo);
00345   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00346     HcalZSThreshold item(cell->rawId(),0);
00347     result->addValues(item);
00348   }
00349   return result;
00350 }
00351 
00352 
00353 std::auto_ptr<HcalL1TriggerObjects> HcalHardcodeCalibrations::produceL1TriggerObjects (const HcalL1TriggerObjectsRcd& rcd) {
00354   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceL1TriggerObjects-> ...";
00355   edm::ESHandle<HcalTopology> htopo;
00356   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00357   const HcalTopology* topo=&(*htopo);
00358 
00359   std::auto_ptr<HcalL1TriggerObjects> result (new HcalL1TriggerObjects (topo));
00360   std::vector <HcalGenericDetId> cells = allCells(*topo);
00361   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00362     HcalL1TriggerObject item(cell->rawId(),0., 1., 0);
00363     result->addValues(item);
00364   }
00365   // add tag and algo values
00366   result->setTagString("hardcoded");
00367   result->setAlgoString("hardcoded");
00368   return result;
00369 }
00370 
00371 
00372 
00373 
00374 std::auto_ptr<HcalElectronicsMap> HcalHardcodeCalibrations::produceElectronicsMap (const HcalElectronicsMapRcd& rcd) {
00375   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceElectronicsMap-> ...";
00376 
00377   std::auto_ptr<HcalElectronicsMap> result (new HcalElectronicsMap ());
00378   HcalDbHardcode::makeHardcodeMap(*result);
00379   return result;
00380 }
00381 
00382 std::auto_ptr<HcalValidationCorrs> HcalHardcodeCalibrations::produceValidationCorrs (const HcalValidationCorrsRcd& rcd) {
00383   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceValidationCorrs-> ...";
00384   edm::ESHandle<HcalTopology> htopo;
00385   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00386   const HcalTopology* topo=&(*htopo);
00387 
00388   std::auto_ptr<HcalValidationCorrs> result (new HcalValidationCorrs (topo));
00389   std::vector <HcalGenericDetId> cells = allCells(*topo);
00390   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00391     HcalValidationCorr item(cell->rawId(),1.0);
00392     result->addValues(item);
00393   }
00394   return result;
00395 }
00396 
00397 std::auto_ptr<HcalLutMetadata> HcalHardcodeCalibrations::produceLutMetadata (const HcalLutMetadataRcd& rcd) {
00398   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLutMetadata-> ...";
00399   edm::ESHandle<HcalTopology> htopo;
00400   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00401   const HcalTopology* topo=&(*htopo);
00402 
00403   std::auto_ptr<HcalLutMetadata> result (new HcalLutMetadata (topo));
00404 
00405   result->setRctLsb( 0.25 );
00406   result->setNominalGain( 0.177 );
00407 
00408   std::vector <HcalGenericDetId> cells = allCells(*topo);
00409   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00410     HcalLutMetadatum item(cell->rawId(),1.0,1,1);
00411     result->addValues(item);
00412   }
00413   return result;
00414 }
00415 
00416 std::auto_ptr<HcalDcsValues> 
00417   HcalHardcodeCalibrations::produceDcsValues (const HcalDcsRcd& rcd) {
00418   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsValues-> ...";
00419   std::auto_ptr<HcalDcsValues> result(new HcalDcsValues);
00420   return result;
00421 }
00422 
00423 std::auto_ptr<HcalDcsMap> HcalHardcodeCalibrations::produceDcsMap (const HcalDcsMapRcd& rcd) {
00424   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsMap-> ...";
00425 
00426   std::auto_ptr<HcalDcsMap> result (new HcalDcsMap ());
00427   HcalDbHardcode::makeHardcodeDcsMap(*result);
00428   return result;
00429 }
00430 
00431 std::auto_ptr<HcalRecoParams> HcalHardcodeCalibrations::produceRecoParams (const HcalRecoParamsRcd& rec) {
00432   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRecoParams-> ...";
00433   edm::ESHandle<HcalTopology> htopo;
00434   rec.getRecord<IdealGeometryRecord>().get(htopo);
00435   const HcalTopology* topo=&(*htopo);
00436 
00437   std::auto_ptr<HcalRecoParams> result (new HcalRecoParams (topo));
00438   std::vector <HcalGenericDetId> cells = allCells(*topo);
00439   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00440     HcalRecoParam item = HcalDbHardcode::makeRecoParam (*cell);
00441     result->addValues(item);
00442   }
00443   return result;
00444 }
00445 std::auto_ptr<HcalTimingParams> HcalHardcodeCalibrations::produceTimingParams (const HcalTimingParamsRcd& rec) {
00446   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimingParams-> ...";
00447   edm::ESHandle<HcalTopology> htopo;
00448   rec.getRecord<IdealGeometryRecord>().get(htopo);
00449   const HcalTopology* topo=&(*htopo);
00450 
00451   std::auto_ptr<HcalTimingParams> result (new HcalTimingParams (topo));
00452   std::vector <HcalGenericDetId> cells = allCells(*topo);
00453   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00454     HcalTimingParam item = HcalDbHardcode::makeTimingParam (*cell);
00455     result->addValues(item);
00456   }
00457   return result;
00458 }
00459 std::auto_ptr<HcalLongRecoParams> HcalHardcodeCalibrations::produceLongRecoParams (const HcalLongRecoParamsRcd& rec) {
00460   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLongRecoParams-> ...";
00461   edm::ESHandle<HcalTopology> htopo;
00462   rec.getRecord<IdealGeometryRecord>().get(htopo);
00463   const HcalTopology* topo=&(*htopo);
00464 
00465   std::auto_ptr<HcalLongRecoParams> result (new HcalLongRecoParams (topo));
00466   std::vector <HcalGenericDetId> cells = allCells(*topo);
00467   std::vector <unsigned int> mSignal; 
00468   mSignal.push_back(4); 
00469   mSignal.push_back(5); 
00470   mSignal.push_back(6);
00471   std::vector <unsigned int> mNoise;  
00472   mNoise.push_back(1);  
00473   mNoise.push_back(2);  
00474   mNoise.push_back(3);
00475   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00476     if (cell->isHcalZDCDetId())
00477       {
00478         HcalLongRecoParam item(cell->rawId(),mSignal,mNoise);
00479         result->addValues(item);
00480       }
00481   }
00482   return result;
00483 }
00484 
00485 std::auto_ptr<HcalMCParams> HcalHardcodeCalibrations::produceMCParams (const HcalMCParamsRcd& rec) {
00486   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceMCParams-> ...";
00487   edm::ESHandle<HcalTopology> htopo;
00488   rec.getRecord<IdealGeometryRecord>().get(htopo);
00489   const HcalTopology* topo=&(*htopo);
00490 
00491   std::auto_ptr<HcalMCParams> result (new HcalMCParams (topo));
00492   std::vector <HcalGenericDetId> cells = allCells(*topo);
00493   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00494     HcalMCParam item(cell->rawId(),0);
00495     result->addValues(item);
00496   }
00497   return result;
00498 }
00499 
00500 
00501 std::auto_ptr<HcalFlagHFDigiTimeParams> HcalHardcodeCalibrations::produceFlagHFDigiTimeParams (const HcalFlagHFDigiTimeParamsRcd& rec) {
00502   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceFlagHFDigiTimeParams-> ...";
00503   edm::ESHandle<HcalTopology> htopo;
00504   rec.getRecord<IdealGeometryRecord>().get(htopo);
00505   const HcalTopology* topo=&(*htopo);
00506 
00507   std::auto_ptr<HcalFlagHFDigiTimeParams> result (new HcalFlagHFDigiTimeParams (topo));
00508   std::vector <HcalGenericDetId> cells = allCells(*topo);
00509   
00510   std::vector<double> coef;
00511   coef.push_back(0.93);
00512   coef.push_back(-0.38275);
00513   coef.push_back(-0.012667);
00514 
00515   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00516     HcalFlagHFDigiTimeParam item(cell->rawId(),
00517                                  1, //firstsample
00518                                  3, // samplestoadd
00519                                  2, //expectedpeak
00520                                  40., // min energy threshold
00521                                  coef // coefficients
00522                                  );
00523     result->addValues(item);
00524   }
00525   return result;
00526 } // produceFlagHFDigiTimeParams;