CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/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.41 2013/04/20 00:36:15 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 
00038   /*
00039   std::cout << std::endl << "HcalHardcodeCalibrations:   maxDepthHB, maxDepthHE = " 
00040             <<  maxDepthHB << ", " <<  maxDepthHE << std::endl;
00041   */
00042 
00043   if (result.size () <= 0) {
00044     for (int eta = -50; eta < 50; eta++) {
00045       for (int phi = 0; phi < 100; phi++) {
00046         for (int depth = 1; depth < maxDepthHB + maxDepthHE; depth++) {
00047           for (int det = 1; det < 5; det++) {
00048             HcalDetId cell ((HcalSubdetector) det, eta, phi, depth);
00049             if (hcaltopology.valid(cell)) result.push_back (cell);
00050 
00051             /*
00052             if (hcaltopology.valid(cell))  
00053               std::cout << " HcalHardcodedCalibrations: det, eta, phi, depth = "
00054                         << det << ",  " << eta << ", " << phi << " , "
00055                         << depth << std::endl;  
00056             */
00057           }
00058         }
00059       }
00060     } 
00061     ZdcTopology zdctopology;
00062     HcalZDCDetId zcell;
00063     HcalZDCDetId::Section section  = HcalZDCDetId::EM;
00064     for(int depth= 1; depth < 6; depth++){
00065       zcell = HcalZDCDetId(section, true, depth);
00066       if(zdctopology.valid(zcell)) result.push_back(zcell);
00067       zcell = HcalZDCDetId(section, false, depth);
00068       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00069      }
00070     section = HcalZDCDetId::HAD;
00071     for(int depth= 1; depth < 5; depth++){
00072       zcell = HcalZDCDetId(section, true, depth);
00073       if(zdctopology.valid(zcell)) result.push_back(zcell);
00074       zcell = HcalZDCDetId(section, false, depth);
00075       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00076     }
00077     section = HcalZDCDetId::LUM;
00078     for(int depth= 1; depth < 3; depth++){
00079       zcell = HcalZDCDetId(section, true, depth);
00080       if(zdctopology.valid(zcell)) result.push_back(zcell);
00081       zcell = HcalZDCDetId(section, false, depth);
00082       if(zdctopology.valid(zcell)) result.push_back(zcell);     
00083     }
00084 
00085     // HcalGenTriggerTower (HcalGenericSubdetector = 5) 
00086     // NASTY HACK !!!
00087     // - As no valid(cell) check found for HcalTrigTowerDetId 
00088     // to create HT cells (ieta=1-28, iphi=1-72)&(ieta=29-32, iphi=1,5,... 69)
00089 
00090     for (int eta = -32; eta <= 32; eta++) {
00091       if(abs(eta) <= 28 && (eta != 0)) {
00092         for (int phi = 1; phi <= 72; phi++) {
00093           HcalTrigTowerDetId cell(eta, phi);       
00094           result.push_back (cell);
00095         }
00096       }
00097       else if (abs(eta) > 28) {
00098         for (int phi = 1; phi <= 69;) {
00099           HcalTrigTowerDetId cell(eta, phi);       
00100           result.push_back (cell);
00101           phi += 4;
00102         }
00103       }
00104     }
00105   }
00106   return result;
00107 }
00108 
00109 }
00110 
00111 HcalHardcodeCalibrations::HcalHardcodeCalibrations ( const edm::ParameterSet& iConfig ): he_recalibration(0), hf_recalibration(0)
00112 {
00113   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::HcalHardcodeCalibrations->...";
00114 
00115   // HE recalibration preparation
00116   iLumi = 0.;
00117   if ( iConfig.exists("iLumi") )
00118     iLumi=iConfig.getParameter<double>("iLumi");
00119 
00120   if( iLumi > 0.0 ) {
00121     bool he_recalib = iConfig.getParameter<bool>("HERecalibration");
00122     bool hf_recalib = iConfig.getParameter<bool>("HFRecalibration");
00123     if(he_recalib)  he_recalibration = new HERecalibration(iLumi);
00124     if(hf_recalib)  hf_recalibration = new HFRecalibration();
00125     
00126   //    std::cout << " HcalHardcodeCalibrations:  iLumi = " <<  iLumi << std::endl;
00127   }
00128 
00129   bool relabel_=false;
00130   edm::ParameterSet ps0;
00131   if ( iConfig.exists("HcalReLabel") ) {
00132     ps0 = iConfig.getParameter<edm::ParameterSet>("HcalReLabel");
00133     relabel_= ps0.getUntrackedParameter<bool>("RelabelHits",false);
00134   }
00135 
00136   if (relabel_) {
00137     std::vector<std::vector<int>> m_segmentation;
00138     m_segmentation.resize(29);
00139     edm::ParameterSet ps1 = ps0.getUntrackedParameter<edm::ParameterSet>("RelabelRules");
00140     for (int i = 0; i < 29; i++) {
00141       char name[10];
00142       snprintf(name,10,"Eta%d",i+1);
00143       if (i > 0) {
00144         m_segmentation[i]=
00145           ps1.getUntrackedParameter<std::vector<int>>(name,m_segmentation[i-1]);
00146       } else {
00147         m_segmentation[i]=ps1.getUntrackedParameter<std::vector<int> >(name);
00148       }
00149       
00150       /*
00151            std::cout << name;
00152       for (unsigned int k=0; k<m_segmentation[i].size(); ++k) {
00153         std::cout << " [" << k << "] " << m_segmentation[i][k];
00154       }
00155       std::cout << std::endl;
00156       */     
00157 
00158     }
00159 
00160     if(he_recalibration !=0) he_recalibration->setDsegm(m_segmentation);
00161   }
00162 
00163 
00164   std::vector <std::string> toGet = iConfig.getUntrackedParameter <std::vector <std::string> > ("toGet");
00165   for(std::vector <std::string>::iterator objectName = toGet.begin(); objectName != toGet.end(); ++objectName ) {
00166     bool all = *objectName == "all";
00167     if ((*objectName == "Pedestals") || all) {
00168       setWhatProduced (this, &HcalHardcodeCalibrations::producePedestals);
00169       findingRecord <HcalPedestalsRcd> ();
00170     }
00171     if ((*objectName == "PedestalWidths") || all) {
00172       setWhatProduced (this, &HcalHardcodeCalibrations::producePedestalWidths);
00173       findingRecord <HcalPedestalWidthsRcd> ();
00174     }
00175     if ((*objectName == "Gains") || all) {
00176       setWhatProduced (this, &HcalHardcodeCalibrations::produceGains);
00177       findingRecord <HcalGainsRcd> ();
00178     }
00179     if ((*objectName == "GainWidths") || all) {
00180       setWhatProduced (this, &HcalHardcodeCalibrations::produceGainWidths);
00181       findingRecord <HcalGainWidthsRcd> ();
00182     }
00183     if ((*objectName == "QIEData") || all) {
00184       setWhatProduced (this, &HcalHardcodeCalibrations::produceQIEData);
00185       findingRecord <HcalQIEDataRcd> ();
00186     }
00187     if ((*objectName == "ChannelQuality") || (*objectName == "channelQuality") || all) {
00188       setWhatProduced (this, &HcalHardcodeCalibrations::produceChannelQuality);
00189       findingRecord <HcalChannelQualityRcd> ();
00190     }
00191     if ((*objectName == "ElectronicsMap") || (*objectName == "electronicsMap") || all) {
00192       setWhatProduced (this, &HcalHardcodeCalibrations::produceElectronicsMap);
00193       findingRecord <HcalElectronicsMapRcd> ();
00194     }
00195     if ((*objectName == "ZSThresholds") || (*objectName == "zsThresholds") || all) {
00196       setWhatProduced (this, &HcalHardcodeCalibrations::produceZSThresholds);
00197       findingRecord <HcalZSThresholdsRcd> ();
00198     }
00199     if ((*objectName == "RespCorrs") || (*objectName == "ResponseCorrection") || all) {
00200       setWhatProduced (this, &HcalHardcodeCalibrations::produceRespCorrs);
00201       findingRecord <HcalRespCorrsRcd> ();
00202     }
00203     if ((*objectName == "LUTCorrs") || (*objectName == "LUTCorrection") || all) {
00204       setWhatProduced (this, &HcalHardcodeCalibrations::produceLUTCorrs);
00205       findingRecord <HcalLUTCorrsRcd> ();
00206     }
00207     if ((*objectName == "PFCorrs") || (*objectName == "PFCorrection") || all) {
00208       setWhatProduced (this, &HcalHardcodeCalibrations::producePFCorrs);
00209       findingRecord <HcalPFCorrsRcd> ();
00210     }
00211     if ((*objectName == "TimeCorrs") || (*objectName == "TimeCorrection") || all) {
00212       setWhatProduced (this, &HcalHardcodeCalibrations::produceTimeCorrs);
00213       findingRecord <HcalTimeCorrsRcd> ();
00214     }
00215     if ((*objectName == "L1TriggerObjects") || (*objectName == "L1Trigger") || all) {
00216       setWhatProduced (this, &HcalHardcodeCalibrations::produceL1TriggerObjects);
00217       findingRecord <HcalL1TriggerObjectsRcd> ();
00218     }
00219     if ((*objectName == "ValidationCorrs") || (*objectName == "ValidationCorrection") || all) {
00220       setWhatProduced (this, &HcalHardcodeCalibrations::produceValidationCorrs);
00221       findingRecord <HcalValidationCorrsRcd> ();
00222     }
00223     if ((*objectName == "LutMetadata") || (*objectName == "lutMetadata") || all) {
00224       setWhatProduced (this, &HcalHardcodeCalibrations::produceLutMetadata);
00225       findingRecord <HcalLutMetadataRcd> ();
00226     }
00227     if ((*objectName == "DcsValues") || all) {
00228       setWhatProduced (this, &HcalHardcodeCalibrations::produceDcsValues);
00229       findingRecord <HcalDcsRcd> ();
00230     }
00231     if ((*objectName == "DcsMap") || (*objectName == "dcsMap") || all) {
00232       setWhatProduced (this, &HcalHardcodeCalibrations::produceDcsMap);
00233       findingRecord <HcalDcsMapRcd> ();
00234     }
00235     if ((*objectName == "RecoParams") || all) {
00236       setWhatProduced (this, &HcalHardcodeCalibrations::produceRecoParams);
00237       findingRecord <HcalRecoParamsRcd> ();
00238     }
00239     if ((*objectName == "LongRecoParams") || all) {
00240       setWhatProduced (this, &HcalHardcodeCalibrations::produceLongRecoParams);
00241       findingRecord <HcalLongRecoParamsRcd> ();
00242     }
00243     if ((*objectName == "MCParams") || all) {
00244       setWhatProduced (this, &HcalHardcodeCalibrations::produceMCParams);
00245       findingRecord <HcalMCParamsRcd> ();
00246     }
00247     if ((*objectName == "FlagHFDigiTimeParams") || all) {
00248       setWhatProduced (this, &HcalHardcodeCalibrations::produceFlagHFDigiTimeParams);
00249       findingRecord <HcalFlagHFDigiTimeParamsRcd> ();
00250     }
00251     if ((*objectName == "CholeskyMatrices") || all) {
00252       setWhatProduced (this, &HcalHardcodeCalibrations::produceCholeskyMatrices);
00253       findingRecord <HcalCholeskyMatricesRcd> ();
00254     }
00255     if ((*objectName == "CovarianceMatrices") || all) {
00256       setWhatProduced (this, &HcalHardcodeCalibrations::produceCovarianceMatrices);
00257       findingRecord <HcalCovarianceMatricesRcd> ();
00258     }
00259   }
00260 }
00261 
00262 
00263 HcalHardcodeCalibrations::~HcalHardcodeCalibrations()
00264 {
00265   if (he_recalibration != 0 ) delete he_recalibration;
00266   if (hf_recalibration != 0 ) delete hf_recalibration;
00267 }
00268 
00269 //
00270 // member functions
00271 //
00272 void 
00273 HcalHardcodeCalibrations::setIntervalFor( const edm::eventsetup::EventSetupRecordKey& iKey, const edm::IOVSyncValue& iTime, edm::ValidityInterval& oInterval ) {
00274   std::string record = iKey.name ();
00275   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::setIntervalFor-> key: " << record << " time: " << iTime.eventID() << '/' << iTime.time ().value ();
00276   oInterval = edm::ValidityInterval (edm::IOVSyncValue::beginOfTime(), edm::IOVSyncValue::endOfTime()); //infinite
00277 }
00278 
00279 std::auto_ptr<HcalPedestals> HcalHardcodeCalibrations::producePedestals (const HcalPedestalsRcd& rec) {
00280   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePedestals-> ...";
00281   edm::ESHandle<HcalTopology> htopo;
00282   rec.getRecord<IdealGeometryRecord>().get(htopo);
00283   const HcalTopology* topo=&(*htopo);
00284 
00285   std::auto_ptr<HcalPedestals> result (new HcalPedestals (topo,false));
00286   std::vector <HcalGenericDetId> cells = allCells(*topo);
00287   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00288     HcalPedestal item = HcalDbHardcode::makePedestal (*cell, false, iLumi);
00289     result->addValues(item);
00290   }
00291   return result;
00292 }
00293 
00294 std::auto_ptr<HcalPedestalWidths> HcalHardcodeCalibrations::producePedestalWidths (const HcalPedestalWidthsRcd& rec) {
00295   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePedestalWidths-> ...";
00296   edm::ESHandle<HcalTopology> htopo;
00297   rec.getRecord<IdealGeometryRecord>().get(htopo);
00298   const HcalTopology* topo=&(*htopo);
00299 
00300   std::auto_ptr<HcalPedestalWidths> result (new HcalPedestalWidths (topo,false));
00301   std::vector <HcalGenericDetId> cells = allCells(*htopo);
00302   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00303     HcalPedestalWidth item = HcalDbHardcode::makePedestalWidth (*cell, iLumi);
00304     result->addValues(item);
00305   }
00306   return result;
00307 }
00308 
00309 std::auto_ptr<HcalGains> HcalHardcodeCalibrations::produceGains (const HcalGainsRcd& rec) {
00310   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGains-> ...";
00311   edm::ESHandle<HcalTopology> htopo;
00312   rec.getRecord<IdealGeometryRecord>().get(htopo);
00313   const HcalTopology* topo=&(*htopo);
00314 
00315   std::auto_ptr<HcalGains> result (new HcalGains (topo));
00316   std::vector <HcalGenericDetId> cells = allCells(*topo);
00317   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00318     HcalGain item = HcalDbHardcode::makeGain (*cell);
00319     result->addValues(item);
00320   }
00321   return result;
00322 }
00323 
00324 std::auto_ptr<HcalGainWidths> HcalHardcodeCalibrations::produceGainWidths (const HcalGainWidthsRcd& rec) {
00325   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGainWidths-> ...";
00326   edm::ESHandle<HcalTopology> htopo;
00327   rec.getRecord<IdealGeometryRecord>().get(htopo);
00328   const HcalTopology* topo=&(*htopo);
00329 
00330   std::auto_ptr<HcalGainWidths> result (new HcalGainWidths (topo));
00331   std::vector <HcalGenericDetId> cells = allCells(*topo);
00332   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00333 
00334     if ( !cell->isHcalTrigTowerDetId()) {
00335       HcalGainWidth item = HcalDbHardcode::makeGainWidth (*cell);
00336       result->addValues(item);
00337     }
00338   }
00339   return result;
00340 }
00341 
00342 std::auto_ptr<HcalQIEData> HcalHardcodeCalibrations::produceQIEData (const HcalQIEDataRcd& rcd) {
00343   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceQIEData-> ...";
00344 
00345   /*
00346   std::cout << std::endl << ">>>  HcalHardcodeCalibrations::produceQIEData"
00347             << std::endl;  
00348   */
00349 
00350   edm::ESHandle<HcalTopology> htopo;
00351   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00352   const HcalTopology* topo=&(*htopo);
00353 
00354   std::auto_ptr<HcalQIEData> result (new HcalQIEData (topo));
00355   std::vector <HcalGenericDetId> cells = allCells(*topo);
00356   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00357     HcalQIECoder coder = HcalDbHardcode::makeQIECoder (*cell);
00358     result->addCoder (coder);
00359   }
00360   return result;
00361 }
00362 
00363 std::auto_ptr<HcalChannelQuality> HcalHardcodeCalibrations::produceChannelQuality (const HcalChannelQualityRcd& rcd) {
00364   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceChannelQuality-> ...";
00365   edm::ESHandle<HcalTopology> htopo;
00366   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00367   const HcalTopology* topo=&(*htopo);
00368 
00369   std::auto_ptr<HcalChannelQuality> result (new HcalChannelQuality (topo));
00370   std::vector <HcalGenericDetId> cells = allCells(*topo);
00371   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00372     HcalChannelStatus item(cell->rawId(),0);
00373     result->addValues(item);
00374   }
00375   return result;
00376 }
00377 
00378 
00379 std::auto_ptr<HcalRespCorrs> HcalHardcodeCalibrations::produceRespCorrs (const HcalRespCorrsRcd& rcd) {
00380   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRespCorrs-> ...";
00381   edm::ESHandle<HcalTopology> htopo;
00382   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00383   const HcalTopology* topo=&(*htopo);
00384 
00385   std::auto_ptr<HcalRespCorrs> result (new HcalRespCorrs (topo));
00386   std::vector <HcalGenericDetId> cells = allCells(*topo);
00387   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00388 
00389     double corr = 1.0; 
00390 
00391     if ((he_recalibration != 0 ) && 
00392         ((*cell).genericSubdet() == HcalGenericDetId::HcalGenEndcap)) {
00393       
00394       int depth_ = HcalDetId(*cell).depth();
00395       int ieta_  = HcalDetId(*cell).ieta();
00396       corr = he_recalibration->getCorr(ieta_, depth_); 
00397       
00398       /*
00399         std::cout << "HE ieta, depth = " << ieta_  << ",  " << depth_  
00400         << "   corr = "  << corr << std::endl;
00401       */
00402 
00403     }
00404     else if ((hf_recalibration != 0 ) && 
00405         ((*cell).genericSubdet() == HcalGenericDetId::HcalGenForward)) {   
00406       int depth_ = HcalDetId(*cell).depth();
00407       int ieta_  = HcalDetId(*cell).ieta();
00408       corr = hf_recalibration->getCorr(ieta_, depth_, iLumi); 
00409 
00410       /*
00411         std::cout << "HF ieta, depth = " << ieta_  << ",  " << depth_  
00412         << "   corr = "  << corr << std::endl;
00413       */
00414 
00415     }
00416 
00417     HcalRespCorr item(cell->rawId(),corr);
00418     result->addValues(item);
00419   }
00420   return result;
00421 }
00422 
00423 std::auto_ptr<HcalLUTCorrs> HcalHardcodeCalibrations::produceLUTCorrs (const HcalLUTCorrsRcd& rcd) {
00424   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLUTCorrs-> ...";
00425   edm::ESHandle<HcalTopology> htopo;
00426   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00427   const HcalTopology* topo=&(*htopo);
00428 
00429   std::auto_ptr<HcalLUTCorrs> result (new HcalLUTCorrs (topo));
00430   std::vector <HcalGenericDetId> cells = allCells(*topo);
00431   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00432     HcalLUTCorr item(cell->rawId(),1.0);
00433     result->addValues(item);
00434   }
00435   return result;
00436 }
00437 
00438 std::auto_ptr<HcalPFCorrs> HcalHardcodeCalibrations::producePFCorrs (const HcalPFCorrsRcd& rcd) {
00439   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePFCorrs-> ...";
00440   edm::ESHandle<HcalTopology> htopo;
00441   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00442   const HcalTopology* topo=&(*htopo);
00443 
00444   std::auto_ptr<HcalPFCorrs> result (new HcalPFCorrs (topo));
00445   std::vector <HcalGenericDetId> cells = allCells(*topo);
00446   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00447     HcalPFCorr item(cell->rawId(),1.0);
00448     result->addValues(item);
00449   }
00450   return result;
00451 }
00452 
00453 std::auto_ptr<HcalTimeCorrs> HcalHardcodeCalibrations::produceTimeCorrs (const HcalTimeCorrsRcd& rcd) {
00454   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimeCorrs-> ...";
00455   edm::ESHandle<HcalTopology> htopo;
00456   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00457   const HcalTopology* topo=&(*htopo);
00458 
00459   std::auto_ptr<HcalTimeCorrs> result (new HcalTimeCorrs (topo));
00460   std::vector <HcalGenericDetId> cells = allCells(*topo);
00461   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00462     HcalTimeCorr item(cell->rawId(),0.0);
00463     result->addValues(item);
00464   }
00465   return result;
00466 }
00467 
00468 std::auto_ptr<HcalZSThresholds> HcalHardcodeCalibrations::produceZSThresholds (const HcalZSThresholdsRcd& rcd) {
00469   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceZSThresholds-> ...";
00470   edm::ESHandle<HcalTopology> htopo;
00471   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00472   const HcalTopology* topo=&(*htopo);
00473 
00474   std::auto_ptr<HcalZSThresholds> result (new HcalZSThresholds (topo));
00475   std::vector <HcalGenericDetId> cells = allCells(*topo);
00476   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00477     HcalZSThreshold item(cell->rawId(),0);
00478     result->addValues(item);
00479   }
00480   return result;
00481 }
00482 
00483 
00484 std::auto_ptr<HcalL1TriggerObjects> HcalHardcodeCalibrations::produceL1TriggerObjects (const HcalL1TriggerObjectsRcd& rcd) {
00485   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceL1TriggerObjects-> ...";
00486   edm::ESHandle<HcalTopology> htopo;
00487   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00488   const HcalTopology* topo=&(*htopo);
00489 
00490   std::auto_ptr<HcalL1TriggerObjects> result (new HcalL1TriggerObjects (topo));
00491   std::vector <HcalGenericDetId> cells = allCells(*topo);
00492   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00493     HcalL1TriggerObject item(cell->rawId(),0., 1., 0);
00494     result->addValues(item);
00495   }
00496   // add tag and algo values
00497   result->setTagString("hardcoded");
00498   result->setAlgoString("hardcoded");
00499   return result;
00500 }
00501 
00502 
00503 
00504 
00505 std::auto_ptr<HcalElectronicsMap> HcalHardcodeCalibrations::produceElectronicsMap (const HcalElectronicsMapRcd& rcd) {
00506   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceElectronicsMap-> ...";
00507 
00508   std::auto_ptr<HcalElectronicsMap> result (new HcalElectronicsMap ());
00509   HcalDbHardcode::makeHardcodeMap(*result);
00510   return result;
00511 }
00512 
00513 std::auto_ptr<HcalValidationCorrs> HcalHardcodeCalibrations::produceValidationCorrs (const HcalValidationCorrsRcd& rcd) {
00514   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceValidationCorrs-> ...";
00515   edm::ESHandle<HcalTopology> htopo;
00516   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00517   const HcalTopology* topo=&(*htopo);
00518 
00519   std::auto_ptr<HcalValidationCorrs> result (new HcalValidationCorrs (topo));
00520   std::vector <HcalGenericDetId> cells = allCells(*topo);
00521   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00522     HcalValidationCorr item(cell->rawId(),1.0);
00523     result->addValues(item);
00524   }
00525   return result;
00526 }
00527 
00528 std::auto_ptr<HcalLutMetadata> HcalHardcodeCalibrations::produceLutMetadata (const HcalLutMetadataRcd& rcd) {
00529   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLutMetadata-> ...";
00530   edm::ESHandle<HcalTopology> htopo;
00531   rcd.getRecord<IdealGeometryRecord>().get(htopo);
00532   const HcalTopology* topo=&(*htopo);
00533 
00534   std::auto_ptr<HcalLutMetadata> result (new HcalLutMetadata (topo));
00535 
00536   result->setRctLsb( 0.5 );
00537   result->setNominalGain(0.003333);  // for HBHE SiPMs
00538 
00539   std::vector <HcalGenericDetId> cells = allCells(*topo);
00540   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00541 
00542     /*
00543     if (cell->isHcalTrigTowerDetId()) {
00544       HcalTrigTowerDetId ht = HcalTrigTowerDetId(*cell);
00545       int ieta = ht.ieta();
00546       int iphi = ht.iphi();
00547       std::cout << " HcalTrigTower cell (ieta,iphi) = " 
00548                <<  ieta << ",  " << iphi << std::endl;
00549     }
00550     */
00551 
00552     HcalLutMetadatum item(cell->rawId(),1.0,1,1);
00553     result->addValues(item);
00554   }
00555   
00556   return result;
00557 }
00558 
00559 std::auto_ptr<HcalDcsValues> 
00560   HcalHardcodeCalibrations::produceDcsValues (const HcalDcsRcd& rcd) {
00561   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsValues-> ...";
00562   std::auto_ptr<HcalDcsValues> result(new HcalDcsValues);
00563   return result;
00564 }
00565 
00566 std::auto_ptr<HcalDcsMap> HcalHardcodeCalibrations::produceDcsMap (const HcalDcsMapRcd& rcd) {
00567   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsMap-> ...";
00568 
00569   std::auto_ptr<HcalDcsMap> result (new HcalDcsMap ());
00570   HcalDbHardcode::makeHardcodeDcsMap(*result);
00571   return result;
00572 }
00573 
00574 std::auto_ptr<HcalRecoParams> HcalHardcodeCalibrations::produceRecoParams (const HcalRecoParamsRcd& rec) {
00575   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRecoParams-> ...";
00576   edm::ESHandle<HcalTopology> htopo;
00577   rec.getRecord<IdealGeometryRecord>().get(htopo);
00578   const HcalTopology* topo=&(*htopo);
00579 
00580   std::auto_ptr<HcalRecoParams> result (new HcalRecoParams (topo));
00581   std::vector <HcalGenericDetId> cells = allCells(*topo);
00582   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00583     HcalRecoParam item = HcalDbHardcode::makeRecoParam (*cell);
00584     result->addValues(item);
00585   }
00586   return result;
00587 }
00588 std::auto_ptr<HcalTimingParams> HcalHardcodeCalibrations::produceTimingParams (const HcalTimingParamsRcd& rec) {
00589   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimingParams-> ...";
00590   edm::ESHandle<HcalTopology> htopo;
00591   rec.getRecord<IdealGeometryRecord>().get(htopo);
00592   const HcalTopology* topo=&(*htopo);
00593 
00594   std::auto_ptr<HcalTimingParams> result (new HcalTimingParams (topo));
00595   std::vector <HcalGenericDetId> cells = allCells(*topo);
00596   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00597     HcalTimingParam item = HcalDbHardcode::makeTimingParam (*cell);
00598     result->addValues(item);
00599   }
00600   return result;
00601 }
00602 std::auto_ptr<HcalLongRecoParams> HcalHardcodeCalibrations::produceLongRecoParams (const HcalLongRecoParamsRcd& rec) {
00603   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLongRecoParams-> ...";
00604   edm::ESHandle<HcalTopology> htopo;
00605   rec.getRecord<IdealGeometryRecord>().get(htopo);
00606   const HcalTopology* topo=&(*htopo);
00607 
00608   std::auto_ptr<HcalLongRecoParams> result (new HcalLongRecoParams (topo));
00609   std::vector <HcalGenericDetId> cells = allCells(*topo);
00610   std::vector <unsigned int> mSignal; 
00611   mSignal.push_back(4); 
00612   mSignal.push_back(5); 
00613   mSignal.push_back(6);
00614   std::vector <unsigned int> mNoise;  
00615   mNoise.push_back(1);  
00616   mNoise.push_back(2);  
00617   mNoise.push_back(3);
00618   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00619     if (cell->isHcalZDCDetId())
00620       {
00621         HcalLongRecoParam item(cell->rawId(),mSignal,mNoise);
00622         result->addValues(item);
00623       }
00624   }
00625   return result;
00626 }
00627 
00628 std::auto_ptr<HcalMCParams> HcalHardcodeCalibrations::produceMCParams (const HcalMCParamsRcd& rec) {
00629 
00630 
00631   //  std::cout << std::endl << " .... HcalHardcodeCalibrations::produceMCParams ->"<< std::endl;
00632 
00633   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceMCParams-> ...";
00634   edm::ESHandle<HcalTopology> htopo;
00635   rec.getRecord<IdealGeometryRecord>().get(htopo);
00636   const HcalTopology* topo=&(*htopo);
00637   std::auto_ptr<HcalMCParams> result (new HcalMCParams (topo));
00638   std::vector <HcalGenericDetId> cells = allCells(*topo);
00639   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00640 
00641     //    HcalMCParam item(cell->rawId(),0);
00642     HcalMCParam item = HcalDbHardcode::makeMCParam (*cell);
00643     result->addValues(item);
00644   }
00645   return result;
00646 }
00647 
00648 
00649 std::auto_ptr<HcalFlagHFDigiTimeParams> HcalHardcodeCalibrations::produceFlagHFDigiTimeParams (const HcalFlagHFDigiTimeParamsRcd& rec) {
00650   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceFlagHFDigiTimeParams-> ...";
00651   edm::ESHandle<HcalTopology> htopo;
00652   rec.getRecord<IdealGeometryRecord>().get(htopo);
00653   const HcalTopology* topo=&(*htopo);
00654 
00655   std::auto_ptr<HcalFlagHFDigiTimeParams> result (new HcalFlagHFDigiTimeParams (topo));
00656   std::vector <HcalGenericDetId> cells = allCells(*topo);
00657   
00658   std::vector<double> coef;
00659   coef.push_back(0.93);
00660   coef.push_back(-0.38275);
00661   coef.push_back(-0.012667);
00662 
00663   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00664     HcalFlagHFDigiTimeParam item(cell->rawId(),
00665                                  1, //firstsample
00666                                  3, // samplestoadd
00667                                  2, //expectedpeak
00668                                  40., // min energy threshold
00669                                  coef // coefficients
00670                                  );
00671     result->addValues(item);
00672   }
00673   return result;
00674 } 
00675 
00676 
00677 std::auto_ptr<HcalCholeskyMatrices> HcalHardcodeCalibrations::produceCholeskyMatrices (const HcalCholeskyMatricesRcd& rec) {
00678 
00679   edm::ESHandle<HcalTopology> htopo;
00680   rec.getRecord<IdealGeometryRecord>().get(htopo);
00681   const HcalTopology* topo=&(*htopo);
00682   std::auto_ptr<HcalCholeskyMatrices> result (new HcalCholeskyMatrices (topo));
00683 
00684   std::vector <HcalGenericDetId> cells = allCells(*topo);
00685   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00686 
00687     int sub = cell->genericSubdet();
00688 
00689     if (sub == HcalGenericDetId::HcalGenBarrel  || 
00690         sub == HcalGenericDetId::HcalGenEndcap  ||
00691         sub == HcalGenericDetId::HcalGenOuter   ||
00692         sub == HcalGenericDetId::HcalGenForward  ) {
00693       HcalCholeskyMatrix item(cell->rawId());
00694       result->addValues(item);
00695     }
00696   }
00697   return result;
00698 
00699 }
00700 std::auto_ptr<HcalCovarianceMatrices> HcalHardcodeCalibrations::produceCovarianceMatrices (const HcalCovarianceMatricesRcd& rec) {
00701 
00702   edm::ESHandle<HcalTopology> htopo;
00703   rec.getRecord<IdealGeometryRecord>().get(htopo);
00704   const HcalTopology* topo=&(*htopo);
00705   std::auto_ptr<HcalCovarianceMatrices> result (new HcalCovarianceMatrices (topo));
00706   std::vector <HcalGenericDetId> cells = allCells(*topo);
00707   for (std::vector <HcalGenericDetId>::const_iterator cell = cells.begin (); cell != cells.end (); cell++) {
00708 
00709     HcalCovarianceMatrix item(cell->rawId());
00710     result->addValues(item);
00711   }
00712   return result;
00713 }