Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 // Original Author:  Fedor Ratnikov
0003 //
0004 //
0005 
0006 #include <memory>
0007 #include <iostream>
0008 #include <string>
0009 #include <vector>
0010 
0011 #include "FWCore/Framework/interface/ValidityInterval.h"
0012 #include "DataFormats/HcalDetId/interface/HcalZDCDetId.h"
0013 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0014 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0015 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0016 #include "DataFormats/HcalDetId/interface/HcalTrigTowerDetId.h"
0017 #include "FWCore/ParameterSet/interface/FileInPath.h"
0018 
0019 #include "CondFormats/DataRecord/interface/HcalAllRcds.h"
0020 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0021 
0022 #include "Geometry/ForwardGeometry/interface/ZdcTopology.h"
0023 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025 
0026 #include "HcalHardcodeCalibrations.h"
0027 
0028 //#define DebugLog
0029 
0030 // class decleration
0031 //
0032 
0033 using namespace cms;
0034 
0035 namespace {
0036 
0037   const std::vector<HcalGenericDetId>& allCells(const HcalTopology& hcaltopology, bool killHE = false) {
0038     static std::vector<HcalGenericDetId> result;
0039     int maxDepth = hcaltopology.maxDepth();
0040 
0041 #ifdef DebugLog
0042     std::cout << std::endl << "HcalHardcodeCalibrations:   maxDepth = " << maxDepth << std::endl;
0043 #endif
0044 
0045     if (result.empty()) {
0046       for (int eta = -HcalDetId::kHcalEtaMask2; eta <= (int)(HcalDetId::kHcalEtaMask2); eta++) {
0047         for (unsigned int phi = 0; phi <= HcalDetId::kHcalPhiMask2; phi++) {
0048           for (int depth = 1; depth <= maxDepth; depth++) {
0049             for (int det = 1; det <= HcalForward; det++) {
0050               HcalDetId cell((HcalSubdetector)det, eta, phi, depth);
0051               if (killHE && HcalEndcap == cell.subdetId())
0052                 continue;
0053               if (hcaltopology.valid(cell)) {
0054                 result.push_back(cell);
0055 #ifdef DebugLog
0056                 std::cout << " HcalHardcodedCalibrations: det|eta|phi|depth = " << det << "|" << eta << "|" << phi
0057                           << "|" << depth << std::endl;
0058 #endif
0059               }
0060             }
0061           }
0062         }
0063       }
0064       ZdcTopology zdctopology;
0065       HcalZDCDetId zcell;
0066       HcalZDCDetId::Section section = HcalZDCDetId::EM;
0067       for (int depth = 1; depth < 6; depth++) {
0068         zcell = HcalZDCDetId(section, true, depth);
0069         if (zdctopology.valid(zcell))
0070           result.push_back(zcell);
0071         zcell = HcalZDCDetId(section, false, depth);
0072         if (zdctopology.valid(zcell))
0073           result.push_back(zcell);
0074       }
0075       section = HcalZDCDetId::HAD;
0076       for (int depth = 1; depth < 5; depth++) {
0077         zcell = HcalZDCDetId(section, true, depth);
0078         if (zdctopology.valid(zcell))
0079           result.push_back(zcell);
0080         zcell = HcalZDCDetId(section, false, depth);
0081         if (zdctopology.valid(zcell))
0082           result.push_back(zcell);
0083       }
0084       section = HcalZDCDetId::LUM;
0085       for (int depth = 1; depth < 3; depth++) {
0086         zcell = HcalZDCDetId(section, true, depth);
0087         if (zdctopology.valid(zcell))
0088           result.push_back(zcell);
0089         zcell = HcalZDCDetId(section, false, depth);
0090         if (zdctopology.valid(zcell))
0091           result.push_back(zcell);
0092       }
0093       section = HcalZDCDetId::RPD;
0094       for (int depth = 1; depth < 17; depth++) {
0095         zcell = HcalZDCDetId(section, true, depth);
0096         if (zdctopology.valid(zcell))
0097           result.push_back(zcell);
0098         zcell = HcalZDCDetId(section, false, depth);
0099         if (zdctopology.valid(zcell))
0100           result.push_back(zcell);
0101       }
0102 
0103       // HcalGenTriggerTower (HcalGenericSubdetector = 5)
0104       // NASTY HACK !!!
0105       // - As no valid(cell) check found for HcalTrigTowerDetId
0106       // to create HT cells (ieta=1-28, iphi=1-72)&(ieta=29-32, iphi=1,5,... 69)
0107 
0108       for (int vers = 0; vers <= HcalTrigTowerDetId::kHcalVersMask; ++vers) {
0109         for (int depth = 0; depth <= HcalTrigTowerDetId::kHcalDepthMask; ++depth) {
0110           for (int eta = -HcalTrigTowerDetId::kHcalEtaMask; eta <= HcalTrigTowerDetId::kHcalEtaMask; eta++) {
0111             for (int phi = 1; phi <= HcalTrigTowerDetId::kHcalPhiMask; phi++) {
0112               HcalTrigTowerDetId cell(eta, phi, depth, vers);
0113               if (hcaltopology.validHT(cell)) {
0114                 result.push_back(cell);
0115 #ifdef DebugLog
0116                 std::cout << " HcalHardcodedCalibrations: eta|phi|depth|vers = " << eta << "|" << phi << "|" << depth
0117                           << "|" << vers << std::endl;
0118 #endif
0119               }
0120             }
0121           }
0122         }
0123       }
0124     }
0125     return result;
0126   }
0127 
0128 }  // namespace
0129 
0130 HcalHardcodeCalibrations::HcalHardcodeCalibrations(const edm::ParameterSet& iConfig)
0131     : hb_recalibration(nullptr),
0132       he_recalibration(nullptr),
0133       hf_recalibration(nullptr),
0134       setHEdsegm(false),
0135       setHBdsegm(false) {
0136   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::HcalHardcodeCalibrations->...";
0137 
0138   if (iConfig.exists("GainWidthsForTrigPrims"))
0139     switchGainWidthsForTrigPrims = iConfig.getParameter<bool>("GainWidthsForTrigPrims");
0140   else
0141     switchGainWidthsForTrigPrims = false;
0142 
0143   //DB helper preparation
0144   dbHardcode.setHB(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("hb")));
0145   dbHardcode.setHE(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("he")));
0146   dbHardcode.setHF(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("hf")));
0147   dbHardcode.setHO(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("ho")));
0148   dbHardcode.setHBUpgrade(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("hbUpgrade")));
0149   dbHardcode.setHEUpgrade(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("heUpgrade")));
0150   dbHardcode.setHFUpgrade(HcalHardcodeParameters(iConfig.getParameter<edm::ParameterSet>("hfUpgrade")));
0151   dbHardcode.useHBUpgrade(iConfig.getParameter<bool>("useHBUpgrade"));
0152   dbHardcode.useHEUpgrade(iConfig.getParameter<bool>("useHEUpgrade"));
0153   dbHardcode.useHFUpgrade(iConfig.getParameter<bool>("useHFUpgrade"));
0154   dbHardcode.useHOUpgrade(iConfig.getParameter<bool>("useHOUpgrade"));
0155   dbHardcode.testHFQIE10(iConfig.getParameter<bool>("testHFQIE10"));
0156   dbHardcode.testHEPlan1(iConfig.getParameter<bool>("testHEPlan1"));
0157   dbHardcode.setKillHE(iConfig.getParameter<bool>("killHE"));
0158   dbHardcode.setSiPMCharacteristics(iConfig.getParameter<std::vector<edm::ParameterSet>>("SiPMCharacteristics"));
0159 
0160   useLayer0Weight = iConfig.getParameter<bool>("useLayer0Weight");
0161   useIeta18depth1 = iConfig.getParameter<bool>("useIeta18depth1");
0162   testHEPlan1 = iConfig.getParameter<bool>("testHEPlan1");
0163   // HB, HE, HF recalibration preparation
0164   iLumi = iConfig.getParameter<double>("iLumi");
0165 
0166   if (iLumi > 0.0) {
0167     bool hb_recalib = iConfig.getParameter<bool>("HBRecalibration");
0168     bool he_recalib = iConfig.getParameter<bool>("HERecalibration");
0169     bool hf_recalib = iConfig.getParameter<bool>("HFRecalibration");
0170     if (hb_recalib) {
0171       hb_recalibration =
0172           std::make_unique<HBHERecalibration>(iLumi,
0173                                               iConfig.getParameter<double>("HBreCalibCutoff"),
0174                                               iConfig.getParameter<edm::FileInPath>("HBmeanenergies").fullPath());
0175     }
0176     if (he_recalib) {
0177       he_recalibration =
0178           std::make_unique<HBHERecalibration>(iLumi,
0179                                               iConfig.getParameter<double>("HEreCalibCutoff"),
0180                                               iConfig.getParameter<edm::FileInPath>("HEmeanenergies").fullPath());
0181     }
0182     if (hf_recalib && !iConfig.getParameter<edm::ParameterSet>("HFRecalParameterBlock").empty())
0183       hf_recalibration =
0184           std::make_unique<HFRecalibration>(iConfig.getParameter<edm::ParameterSet>("HFRecalParameterBlock"));
0185 
0186 #ifdef DebugLog
0187     std::cout << " HcalHardcodeCalibrations:  iLumi = " << iLumi << std::endl;
0188 #endif
0189   }
0190 
0191   std::vector<std::string> toGet = iConfig.getUntrackedParameter<std::vector<std::string>>("toGet");
0192   for (auto& objectName : toGet) {
0193     bool all = objectName == "all";
0194 #ifdef DebugLog
0195     std::cout << "Load parameters for " << objectName << std::endl;
0196 #endif
0197     if ((objectName == "Pedestals") || all) {
0198       topoTokens_[kPedestals] = setWhatProduced(this, &HcalHardcodeCalibrations::producePedestals).consumes();
0199       findingRecord<HcalPedestalsRcd>();
0200     }
0201     if ((objectName == "PedestalWidths") || all) {
0202       topoTokens_[kPedestalWidths] = setWhatProduced(this, &HcalHardcodeCalibrations::producePedestalWidths).consumes();
0203       findingRecord<HcalPedestalWidthsRcd>();
0204     }
0205     if ((objectName == "EffectivePedestals") || all) {
0206       topoTokens_[kEffectivePedestals] =
0207           setWhatProduced(this, &HcalHardcodeCalibrations::produceEffectivePedestals, edm::es::Label("effective"))
0208               .consumes();
0209       findingRecord<HcalPedestalsRcd>();
0210     }
0211     if ((objectName == "EffectivePedestalWidths") || all) {
0212       topoTokens_[kEffectivePedestalWidths] =
0213           setWhatProduced(this, &HcalHardcodeCalibrations::produceEffectivePedestalWidths, edm::es::Label("effective"))
0214               .consumes();
0215       findingRecord<HcalPedestalWidthsRcd>();
0216     }
0217     if ((objectName == "Gains") || all) {
0218       topoTokens_[kGains] = setWhatProduced(this, &HcalHardcodeCalibrations::produceGains).consumes();
0219       findingRecord<HcalGainsRcd>();
0220     }
0221     if ((objectName == "GainWidths") || all) {
0222       topoTokens_[kGainWidths] = setWhatProduced(this, &HcalHardcodeCalibrations::produceGainWidths).consumes();
0223       findingRecord<HcalGainWidthsRcd>();
0224     }
0225     if ((objectName == "PFCuts") || all) {
0226       topoTokens_[kPFCuts] = setWhatProduced(this, &HcalHardcodeCalibrations::producePFCuts).consumes();
0227       findingRecord<HcalPFCutsRcd>();
0228     }
0229     if ((objectName == "QIEData") || all) {
0230       topoTokens_[kQIEData] = setWhatProduced(this, &HcalHardcodeCalibrations::produceQIEData).consumes();
0231       findingRecord<HcalQIEDataRcd>();
0232     }
0233     if ((objectName == "QIETypes") || all) {
0234       topoTokens_[kQIETypes] = setWhatProduced(this, &HcalHardcodeCalibrations::produceQIETypes).consumes();
0235       findingRecord<HcalQIETypesRcd>();
0236     }
0237     if ((objectName == "ChannelQuality") || (objectName == "channelQuality") || all) {
0238       topoTokens_[kChannelQuality] = setWhatProduced(this, &HcalHardcodeCalibrations::produceChannelQuality).consumes();
0239       findingRecord<HcalChannelQualityRcd>();
0240     }
0241     if ((objectName == "ElectronicsMap") || (objectName == "electronicsMap") || all) {
0242       topoTokens_[kElectronicsMap] = setWhatProduced(this, &HcalHardcodeCalibrations::produceElectronicsMap).consumes();
0243       findingRecord<HcalElectronicsMapRcd>();
0244     }
0245     if ((objectName == "ZSThresholds") || (objectName == "zsThresholds") || all) {
0246       topoTokens_[kZSThresholds] = setWhatProduced(this, &HcalHardcodeCalibrations::produceZSThresholds).consumes();
0247       findingRecord<HcalZSThresholdsRcd>();
0248     }
0249     if ((objectName == "RespCorrs") || (objectName == "ResponseCorrection") || all) {
0250       auto c = setWhatProduced(this, &HcalHardcodeCalibrations::produceRespCorrs);
0251       topoTokens_[kRespCorrs] = c.consumes();
0252       if (he_recalibration) {
0253         heDarkeningToken_ = c.consumes(edm::ESInputTag("", "HE"));
0254       }
0255       if (hb_recalibration) {
0256         hbDarkeningToken_ = c.consumes(edm::ESInputTag("", "HB"));
0257       }
0258       findingRecord<HcalRespCorrsRcd>();
0259     }
0260     if ((objectName == "LUTCorrs") || (objectName == "LUTCorrection") || all) {
0261       topoTokens_[kLUTCorrs] = setWhatProduced(this, &HcalHardcodeCalibrations::produceLUTCorrs).consumes();
0262       findingRecord<HcalLUTCorrsRcd>();
0263     }
0264     if ((objectName == "PFCorrs") || (objectName == "PFCorrection") || all) {
0265       topoTokens_[kPFCorrs] = setWhatProduced(this, &HcalHardcodeCalibrations::producePFCorrs).consumes();
0266       findingRecord<HcalPFCorrsRcd>();
0267     }
0268     if ((objectName == "TimeCorrs") || (objectName == "TimeCorrection") || all) {
0269       topoTokens_[kTimeCorrs] = setWhatProduced(this, &HcalHardcodeCalibrations::produceTimeCorrs).consumes();
0270       findingRecord<HcalTimeCorrsRcd>();
0271     }
0272     if ((objectName == "L1TriggerObjects") || (objectName == "L1Trigger") || all) {
0273       topoTokens_[kL1TriggerObjects] =
0274           setWhatProduced(this, &HcalHardcodeCalibrations::produceL1TriggerObjects).consumes();
0275       findingRecord<HcalL1TriggerObjectsRcd>();
0276     }
0277     if ((objectName == "ValidationCorrs") || (objectName == "ValidationCorrection") || all) {
0278       topoTokens_[kValidationCorrs] =
0279           setWhatProduced(this, &HcalHardcodeCalibrations::produceValidationCorrs).consumes();
0280       findingRecord<HcalValidationCorrsRcd>();
0281     }
0282     if ((objectName == "LutMetadata") || (objectName == "lutMetadata") || all) {
0283       topoTokens_[kLutMetadata] = setWhatProduced(this, &HcalHardcodeCalibrations::produceLutMetadata).consumes();
0284       findingRecord<HcalLutMetadataRcd>();
0285     }
0286     if ((objectName == "DcsValues") || all) {
0287       setWhatProduced(this, &HcalHardcodeCalibrations::produceDcsValues);
0288       findingRecord<HcalDcsRcd>();
0289     }
0290     if ((objectName == "DcsMap") || (objectName == "dcsMap") || all) {
0291       setWhatProduced(this, &HcalHardcodeCalibrations::produceDcsMap);
0292       findingRecord<HcalDcsMapRcd>();
0293     }
0294     if ((objectName == "RecoParams") || all) {
0295       topoTokens_[kRecoParams] = setWhatProduced(this, &HcalHardcodeCalibrations::produceRecoParams).consumes();
0296       findingRecord<HcalRecoParamsRcd>();
0297     }
0298     if ((objectName == "LongRecoParams") || all) {
0299       topoTokens_[kLongRecoParams] = setWhatProduced(this, &HcalHardcodeCalibrations::produceLongRecoParams).consumes();
0300       findingRecord<HcalLongRecoParamsRcd>();
0301     }
0302     if ((objectName == "ZDCLowGainFractions") || all) {
0303       topoTokens_[kZDCLowGainFractions] =
0304           setWhatProduced(this, &HcalHardcodeCalibrations::produceZDCLowGainFractions).consumes();
0305       findingRecord<HcalZDCLowGainFractionsRcd>();
0306     }
0307     if ((objectName == "MCParams") || all) {
0308       topoTokens_[kMCParams] = setWhatProduced(this, &HcalHardcodeCalibrations::produceMCParams).consumes();
0309       findingRecord<HcalMCParamsRcd>();
0310     }
0311     if ((objectName == "FlagHFDigiTimeParams") || all) {
0312       topoTokens_[kFlagHFDigiTimeParams] =
0313           setWhatProduced(this, &HcalHardcodeCalibrations::produceFlagHFDigiTimeParams).consumes();
0314       findingRecord<HcalFlagHFDigiTimeParamsRcd>();
0315     }
0316     if ((objectName == "FrontEndMap") || (objectName == "frontEndMap") || all) {
0317       topoTokens_[kFrontEndMap] = setWhatProduced(this, &HcalHardcodeCalibrations::produceFrontEndMap).consumes();
0318       findingRecord<HcalFrontEndMapRcd>();
0319     }
0320     if ((objectName == "SiPMParameters") || all) {
0321       topoTokens_[kSiPMParameters] = setWhatProduced(this, &HcalHardcodeCalibrations::produceSiPMParameters).consumes();
0322       findingRecord<HcalSiPMParametersRcd>();
0323     }
0324     if ((objectName == "SiPMCharacteristics") || all) {
0325       setWhatProduced(this, &HcalHardcodeCalibrations::produceSiPMCharacteristics);
0326       findingRecord<HcalSiPMCharacteristicsRcd>();
0327     }
0328     if ((objectName == "TPChannelParameters") || all) {
0329       topoTokens_[kTPChannelParameters] =
0330           setWhatProduced(this, &HcalHardcodeCalibrations::produceTPChannelParameters).consumes();
0331       findingRecord<HcalTPChannelParametersRcd>();
0332     }
0333     if ((objectName == "TPParameters") || all) {
0334       setWhatProduced(this, &HcalHardcodeCalibrations::produceTPParameters);
0335       findingRecord<HcalTPParametersRcd>();
0336     }
0337   }
0338 }
0339 
0340 HcalHardcodeCalibrations::~HcalHardcodeCalibrations() {}
0341 
0342 //
0343 // member functions
0344 //
0345 void HcalHardcodeCalibrations::setIntervalFor(const edm::eventsetup::EventSetupRecordKey& iKey,
0346                                               const edm::IOVSyncValue& iTime,
0347                                               edm::ValidityInterval& oInterval) {
0348   std::string record = iKey.name();
0349   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::setIntervalFor-> key: " << record << " time: " << iTime.eventID()
0350                        << '/' << iTime.time().value();
0351   oInterval = edm::ValidityInterval(edm::IOVSyncValue::beginOfTime(), edm::IOVSyncValue::endOfTime());  //infinite
0352 }
0353 
0354 std::unique_ptr<HcalPedestals> HcalHardcodeCalibrations::producePedestals_(
0355     const HcalPedestalsRcd& rec, const edm::ESGetToken<HcalTopology, HcalRecNumberingRecord>& token, bool eff) {
0356   std::string seff = eff ? "Effective" : "";
0357   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produce" << seff << "Pedestals-> ...";
0358 
0359   auto const& topo = rec.get(token);
0360   auto result = std::make_unique<HcalPedestals>(&topo, false);
0361   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0362   for (auto cell : cells) {
0363     HcalPedestal item = dbHardcode.makePedestal(cell, false, eff, &topo, iLumi);
0364     result->addValues(item);
0365   }
0366   return result;
0367 }
0368 
0369 std::unique_ptr<HcalPedestalWidths> HcalHardcodeCalibrations::producePedestalWidths_(
0370     const HcalPedestalWidthsRcd& rec, const edm::ESGetToken<HcalTopology, HcalRecNumberingRecord>& token, bool eff) {
0371   std::string seff = eff ? "Effective" : "";
0372   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produce" << seff << "PedestalWidths-> ...";
0373   auto const& topo = rec.get(token);
0374   auto result = std::make_unique<HcalPedestalWidths>(&topo, false);
0375   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0376   for (auto cell : cells) {
0377     HcalPedestalWidth item = dbHardcode.makePedestalWidth(cell, eff, &topo, iLumi);
0378     result->addValues(item);
0379   }
0380   return result;
0381 }
0382 
0383 std::unique_ptr<HcalPedestals> HcalHardcodeCalibrations::producePedestals(const HcalPedestalsRcd& rec) {
0384   return producePedestals_(rec, topoTokens_[kPedestals], false);
0385 }
0386 
0387 std::unique_ptr<HcalPedestals> HcalHardcodeCalibrations::produceEffectivePedestals(const HcalPedestalsRcd& rec) {
0388   return producePedestals_(rec, topoTokens_[kEffectivePedestals], true);
0389 }
0390 
0391 std::unique_ptr<HcalPedestalWidths> HcalHardcodeCalibrations::producePedestalWidths(const HcalPedestalWidthsRcd& rec) {
0392   return producePedestalWidths_(rec, topoTokens_[kPedestalWidths], false);
0393 }
0394 
0395 std::unique_ptr<HcalPedestalWidths> HcalHardcodeCalibrations::produceEffectivePedestalWidths(
0396     const HcalPedestalWidthsRcd& rec) {
0397   return producePedestalWidths_(rec, topoTokens_[kEffectivePedestalWidths], true);
0398 }
0399 
0400 std::unique_ptr<HcalGains> HcalHardcodeCalibrations::produceGains(const HcalGainsRcd& rec) {
0401   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGains-> ...";
0402 
0403   auto const& topo = rec.get(topoTokens_[kGains]);
0404   auto result = std::make_unique<HcalGains>(&topo);
0405   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0406   for (auto cell : cells) {
0407     HcalGain item = dbHardcode.makeGain(cell);
0408     result->addValues(item);
0409   }
0410   return result;
0411 }
0412 
0413 std::unique_ptr<HcalGainWidths> HcalHardcodeCalibrations::produceGainWidths(const HcalGainWidthsRcd& rec) {
0414   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceGainWidths-> ...";
0415 
0416   auto const& topo = rec.get(topoTokens_[kGainWidths]);
0417   auto result = std::make_unique<HcalGainWidths>(&topo);
0418   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0419   for (auto cell : cells) {
0420     // for Upgrade - include TrigPrims, for regular case - only HcalDetId
0421     if (switchGainWidthsForTrigPrims) {
0422       HcalGainWidth item = dbHardcode.makeGainWidth(cell);
0423       result->addValues(item);
0424     } else if (!cell.isHcalTrigTowerDetId()) {
0425       HcalGainWidth item = dbHardcode.makeGainWidth(cell);
0426       result->addValues(item);
0427     }
0428   }
0429   return result;
0430 }
0431 
0432 std::unique_ptr<HcalPFCuts> HcalHardcodeCalibrations::producePFCuts(const HcalPFCutsRcd& rec) {
0433   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePFCuts-> ...";
0434 
0435   auto const& topo = rec.get(topoTokens_[kPFCuts]);
0436   auto result = std::make_unique<HcalPFCuts>(&topo);
0437   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0438   for (auto cell : cells) {
0439     // Use only standard Hcal channels for now, no TrigPrims
0440     if (!cell.isHcalTrigTowerDetId()) {
0441       HcalPFCut item = dbHardcode.makePFCut(cell, iLumi, dbHardcode.killHE());
0442       result->addValues(item);
0443     }
0444   }
0445   return result;
0446 }
0447 
0448 std::unique_ptr<HcalQIEData> HcalHardcodeCalibrations::produceQIEData(const HcalQIEDataRcd& rcd) {
0449   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceQIEData-> ...";
0450 
0451   /*
0452   std::cout << std::endl << ">>>  HcalHardcodeCalibrations::produceQIEData"
0453         << std::endl;  
0454   */
0455 
0456   auto const& topo = rcd.get(topoTokens_[kQIEData]);
0457   auto result = std::make_unique<HcalQIEData>(&topo);
0458   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0459   for (auto cell : cells) {
0460     HcalQIECoder coder = dbHardcode.makeQIECoder(cell);
0461     result->addCoder(coder);
0462   }
0463   return result;
0464 }
0465 
0466 std::unique_ptr<HcalQIETypes> HcalHardcodeCalibrations::produceQIETypes(const HcalQIETypesRcd& rcd) {
0467   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceQIETypes-> ...";
0468   auto const& topo = rcd.get(topoTokens_[kQIETypes]);
0469 
0470   auto result = std::make_unique<HcalQIETypes>(&topo);
0471   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0472   for (auto cell : cells) {
0473     HcalQIEType item = dbHardcode.makeQIEType(cell);
0474     result->addValues(item);
0475   }
0476   return result;
0477 }
0478 
0479 std::unique_ptr<HcalChannelQuality> HcalHardcodeCalibrations::produceChannelQuality(const HcalChannelQualityRcd& rcd) {
0480   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceChannelQuality-> ...";
0481   auto const& topo = rcd.get(topoTokens_[kChannelQuality]);
0482 
0483   auto result = std::make_unique<HcalChannelQuality>(&topo);
0484   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0485   for (auto cell : cells) {
0486     // Special: removal of (non-instrumented) layer "-1"("nose") = depth 1
0487     // from Upgrade HE, either from
0488     // (i)  HEP17 sector in 2017 or
0489     // (ii) the entire HE rin=18 from 2018 through Run 3.
0490     // May require a revision  by 2021.
0491 
0492     uint32_t status = 0;
0493 
0494     if (!(cell.isHcalZDCDetId())) {
0495       HcalDetId hid = HcalDetId(cell);
0496       int iphi = hid.iphi();
0497       int ieta = hid.ieta();
0498       int absieta = hid.ietaAbs();
0499       int depth = hid.depth();
0500 
0501       // specific HEP17 sector (2017 only)
0502       bool isHEP17 = (iphi >= 63) && (iphi <= 66) && (ieta > 0);
0503       // |ieta|=18, depth=1
0504       bool is18d1 = (absieta == 18) && (depth == 1);
0505 
0506       if ((!useIeta18depth1 && is18d1) && ((testHEPlan1 && isHEP17) || (!testHEPlan1))) {
0507         status = 0x8002;  // dead cell
0508       }
0509     }
0510 
0511     HcalChannelStatus item(cell.rawId(), status);
0512     result->addValues(item);
0513   }
0514 
0515   return result;
0516 }
0517 
0518 std::unique_ptr<HcalRespCorrs> HcalHardcodeCalibrations::produceRespCorrs(const HcalRespCorrsRcd& rcd) {
0519   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRespCorrs-> ...";
0520   auto const& topo = rcd.get(topoTokens_[kRespCorrs]);
0521 
0522   //set depth segmentation for HB/HE recalib - only happens once
0523   if ((he_recalibration && !setHEdsegm) || (hb_recalibration && !setHBdsegm)) {
0524     std::vector<std::vector<int>> m_segmentation;
0525     int maxEta = topo.lastHBHERing();
0526     m_segmentation.resize(maxEta);
0527     for (int i = 0; i < maxEta; i++) {
0528       topo.getDepthSegmentation(i + 1, m_segmentation[i]);
0529     }
0530     if (he_recalibration && !setHEdsegm) {
0531       he_recalibration->setup(m_segmentation, &rcd.get(heDarkeningToken_));
0532       setHEdsegm = true;
0533     }
0534     if (hb_recalibration && !setHBdsegm) {
0535       hb_recalibration->setup(m_segmentation, &rcd.get(hbDarkeningToken_));
0536       setHBdsegm = true;
0537     }
0538   }
0539 
0540   auto result = std::make_unique<HcalRespCorrs>(&topo);
0541   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0542   for (const auto& cell : cells) {
0543     double corr = 1.0;
0544 
0545     //check for layer 0 reweighting: when depth 1 has only one layer, it is layer 0
0546     if (useLayer0Weight &&
0547         ((cell.genericSubdet() == HcalGenericDetId::HcalGenEndcap) ||
0548          (cell.genericSubdet() == HcalGenericDetId::HcalGenBarrel)) &&
0549         (HcalDetId(cell).depth() == 1 &&
0550          dbHardcode.getLayersInDepth(HcalDetId(cell).ietaAbs(), HcalDetId(cell).depth(), &topo) == 1)) {
0551       //layer 0 is thicker than other layers (9mm vs 3.7mm) and brighter (Bicron vs SCSN81)
0552       //in Run1/Run2 (pre-2017 for HE), ODU for layer 0 had neutral density filter attached
0553       //NDF was simulated as weight of 0.5 applied to Geant energy deposits
0554       //for Phase1, NDF is removed - simulated as weight of 1.2 applied to Geant energy deposits
0555       //to maintain RECO calibrations, move the layer 0 energy scale back to its previous state using respcorrs
0556       corr = 0.5 / 1.2;
0557     }
0558 
0559     if ((hb_recalibration != nullptr) && (cell.genericSubdet() == HcalGenericDetId::HcalGenBarrel)) {
0560       int depth_ = HcalDetId(cell).depth();
0561       int ieta_ = HcalDetId(cell).ieta();
0562       corr *= hb_recalibration->getCorr(ieta_, depth_);
0563 #ifdef DebugLog
0564       std::cout << "HB ieta, depth = " << ieta_ << ",  " << depth_ << "   corr = " << corr << std::endl;
0565 #endif
0566     } else if ((he_recalibration != nullptr) && (cell.genericSubdet() == HcalGenericDetId::HcalGenEndcap)) {
0567       int depth_ = HcalDetId(cell).depth();
0568       int ieta_ = HcalDetId(cell).ieta();
0569       corr *= he_recalibration->getCorr(ieta_, depth_);
0570 #ifdef DebugLog
0571       std::cout << "HE ieta, depth = " << ieta_ << ",  " << depth_ << "   corr = " << corr << std::endl;
0572 #endif
0573     } else if ((hf_recalibration != nullptr) && (cell.genericSubdet() == HcalGenericDetId::HcalGenForward)) {
0574       int depth_ = HcalDetId(cell).depth();
0575       int ieta_ = HcalDetId(cell).ieta();
0576       corr = hf_recalibration->getCorr(ieta_, depth_, iLumi);
0577 #ifdef DebugLog
0578       std::cout << "HF ieta, depth = " << ieta_ << ",  " << depth_ << "   corr = " << corr << std::endl;
0579 #endif
0580     }
0581 
0582     HcalRespCorr item(cell.rawId(), corr);
0583     result->addValues(item);
0584   }
0585   return result;
0586 }
0587 
0588 std::unique_ptr<HcalLUTCorrs> HcalHardcodeCalibrations::produceLUTCorrs(const HcalLUTCorrsRcd& rcd) {
0589   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLUTCorrs-> ...";
0590   auto const& topo = rcd.get(topoTokens_[kLUTCorrs]);
0591 
0592   auto result = std::make_unique<HcalLUTCorrs>(&topo);
0593   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0594   for (auto cell : cells) {
0595     HcalLUTCorr item(cell.rawId(), 1.0);
0596     result->addValues(item);
0597   }
0598   return result;
0599 }
0600 
0601 std::unique_ptr<HcalPFCorrs> HcalHardcodeCalibrations::producePFCorrs(const HcalPFCorrsRcd& rcd) {
0602   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::producePFCorrs-> ...";
0603   auto const& topo = rcd.get(topoTokens_[kPFCorrs]);
0604 
0605   auto result = std::make_unique<HcalPFCorrs>(&topo);
0606   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0607   for (auto cell : cells) {
0608     HcalPFCorr item(cell.rawId(), 1.0);
0609     result->addValues(item);
0610   }
0611   return result;
0612 }
0613 
0614 std::unique_ptr<HcalTimeCorrs> HcalHardcodeCalibrations::produceTimeCorrs(const HcalTimeCorrsRcd& rcd) {
0615   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimeCorrs-> ...";
0616   auto const& topo = rcd.get(topoTokens_[kTimeCorrs]);
0617 
0618   auto result = std::make_unique<HcalTimeCorrs>(&topo);
0619   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0620   for (auto cell : cells) {
0621     HcalTimeCorr item(cell.rawId(), 0.0);
0622     result->addValues(item);
0623   }
0624   return result;
0625 }
0626 
0627 std::unique_ptr<HcalZSThresholds> HcalHardcodeCalibrations::produceZSThresholds(const HcalZSThresholdsRcd& rcd) {
0628   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceZSThresholds-> ...";
0629   auto const& topo = rcd.get(topoTokens_[kZSThresholds]);
0630 
0631   auto result = std::make_unique<HcalZSThresholds>(&topo);
0632   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0633   for (auto cell : cells) {
0634     HcalZSThreshold item = dbHardcode.makeZSThreshold(cell);
0635     result->addValues(item);
0636   }
0637   return result;
0638 }
0639 
0640 std::unique_ptr<HcalL1TriggerObjects> HcalHardcodeCalibrations::produceL1TriggerObjects(
0641     const HcalL1TriggerObjectsRcd& rcd) {
0642   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceL1TriggerObjects-> ...";
0643   auto const& topo = rcd.get(topoTokens_[kL1TriggerObjects]);
0644 
0645   auto result = std::make_unique<HcalL1TriggerObjects>(&topo);
0646   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0647   for (auto cell : cells) {
0648     HcalL1TriggerObject item(cell.rawId(), 0., 1., 0);
0649     result->addValues(item);
0650   }
0651   // add tag and algo values
0652   result->setTagString("hardcoded");
0653   result->setAlgoString("hardcoded");
0654   return result;
0655 }
0656 
0657 std::unique_ptr<HcalElectronicsMap> HcalHardcodeCalibrations::produceElectronicsMap(const HcalElectronicsMapRcd& rcd) {
0658   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceElectronicsMap-> ...";
0659   auto const& topo = rcd.get(topoTokens_[kElectronicsMap]);
0660 
0661   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0662   return dbHardcode.makeHardcodeMap(cells);
0663 }
0664 
0665 std::unique_ptr<HcalValidationCorrs> HcalHardcodeCalibrations::produceValidationCorrs(
0666     const HcalValidationCorrsRcd& rcd) {
0667   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceValidationCorrs-> ...";
0668   auto const& topo = rcd.get(topoTokens_[kValidationCorrs]);
0669 
0670   auto result = std::make_unique<HcalValidationCorrs>(&topo);
0671   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0672   for (auto cell : cells) {
0673     HcalValidationCorr item(cell.rawId(), 1.0);
0674     result->addValues(item);
0675   }
0676   return result;
0677 }
0678 
0679 std::unique_ptr<HcalLutMetadata> HcalHardcodeCalibrations::produceLutMetadata(const HcalLutMetadataRcd& rcd) {
0680   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLutMetadata-> ...";
0681   auto const& topo = rcd.get(topoTokens_[kLutMetadata]);
0682 
0683   auto result = std::make_unique<HcalLutMetadata>(&topo);
0684 
0685   result->setRctLsb(0.5);
0686   result->setNominalGain(0.177);  // for HBHE SiPMs
0687 
0688   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0689   for (const auto& cell : cells) {
0690     float rcalib = 1.;
0691     int granularity = 1;
0692     int threshold = 0;
0693 
0694     if (cell.isHcalTrigTowerDetId()) {
0695       rcalib = 0.;
0696     }
0697 
0698     HcalLutMetadatum item(cell.rawId(), rcalib, granularity, threshold);
0699     result->addValues(item);
0700   }
0701 
0702   return result;
0703 }
0704 
0705 std::unique_ptr<HcalDcsValues> HcalHardcodeCalibrations::produceDcsValues(const HcalDcsRcd& rcd) {
0706   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsValues-> ...";
0707   auto result = std::make_unique<HcalDcsValues>();
0708   return result;
0709 }
0710 
0711 std::unique_ptr<HcalDcsMap> HcalHardcodeCalibrations::produceDcsMap(const HcalDcsMapRcd& rcd) {
0712   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceDcsMap-> ...";
0713 
0714   return dbHardcode.makeHardcodeDcsMap();
0715 }
0716 
0717 std::unique_ptr<HcalRecoParams> HcalHardcodeCalibrations::produceRecoParams(const HcalRecoParamsRcd& rec) {
0718   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceRecoParams-> ...";
0719   auto const& topo = rec.get(topoTokens_[kRecoParams]);
0720 
0721   auto result = std::make_unique<HcalRecoParams>(&topo);
0722   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0723   for (auto cell : cells) {
0724     HcalRecoParam item = dbHardcode.makeRecoParam(cell);
0725     result->addValues(item);
0726   }
0727   return result;
0728 }
0729 
0730 std::unique_ptr<HcalTimingParams> HcalHardcodeCalibrations::produceTimingParams(const HcalTimingParamsRcd& rec) {
0731   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTimingParams-> ...";
0732   auto const& topo = rec.get(topoTokens_[kTimingParams]);
0733 
0734   auto result = std::make_unique<HcalTimingParams>(&topo);
0735   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0736   for (auto cell : cells) {
0737     HcalTimingParam item = dbHardcode.makeTimingParam(cell);
0738     result->addValues(item);
0739   }
0740   return result;
0741 }
0742 
0743 std::unique_ptr<HcalLongRecoParams> HcalHardcodeCalibrations::produceLongRecoParams(const HcalLongRecoParamsRcd& rec) {
0744   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceLongRecoParams-> ...";
0745   auto const& topo = rec.get(topoTokens_[kLongRecoParams]);
0746 
0747   auto result = std::make_unique<HcalLongRecoParams>(&topo);
0748   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0749   std::vector<unsigned int> mSignal;
0750   mSignal.push_back(4);
0751   mSignal.push_back(5);
0752   mSignal.push_back(6);
0753   std::vector<unsigned int> mNoise;
0754   mNoise.push_back(1);
0755   mNoise.push_back(2);
0756   mNoise.push_back(3);
0757   for (auto cell : cells) {
0758     if (cell.isHcalZDCDetId()) {
0759       HcalLongRecoParam item(cell.rawId(), mSignal, mNoise);
0760       result->addValues(item);
0761     }
0762   }
0763   return result;
0764 }
0765 
0766 std::unique_ptr<HcalZDCLowGainFractions> HcalHardcodeCalibrations::produceZDCLowGainFractions(
0767     const HcalZDCLowGainFractionsRcd& rec) {
0768   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceZDCLowGainFractions-> ...";
0769   auto const& topo = rec.get(topoTokens_[kZDCLowGainFractions]);
0770 
0771   auto result = std::make_unique<HcalZDCLowGainFractions>(&topo);
0772   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0773   for (auto cell : cells) {
0774     HcalZDCLowGainFraction item(cell.rawId(), 0.0);
0775     result->addValues(item);
0776   }
0777   return result;
0778 }
0779 
0780 std::unique_ptr<HcalMCParams> HcalHardcodeCalibrations::produceMCParams(const HcalMCParamsRcd& rec) {
0781   //  std::cout << std::endl << " .... HcalHardcodeCalibrations::produceMCParams ->"<< std::endl;
0782 
0783   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceMCParams-> ...";
0784   auto const& topo = rec.get(topoTokens_[kMCParams]);
0785   auto result = std::make_unique<HcalMCParams>(&topo);
0786   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0787   for (auto cell : cells) {
0788     HcalMCParam item = dbHardcode.makeMCParam(cell);
0789     result->addValues(item);
0790   }
0791   return result;
0792 }
0793 
0794 std::unique_ptr<HcalFlagHFDigiTimeParams> HcalHardcodeCalibrations::produceFlagHFDigiTimeParams(
0795     const HcalFlagHFDigiTimeParamsRcd& rec) {
0796   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceFlagHFDigiTimeParams-> ...";
0797   auto const& topo = rec.get(topoTokens_[kFlagHFDigiTimeParams]);
0798 
0799   auto result = std::make_unique<HcalFlagHFDigiTimeParams>(&topo);
0800   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0801 
0802   std::vector<double> coef;
0803   coef.push_back(0.93);
0804   coef.push_back(-0.38275);
0805   coef.push_back(-0.012667);
0806 
0807   for (auto cell : cells) {
0808     HcalFlagHFDigiTimeParam item(cell.rawId(),
0809                                  1,    //firstsample
0810                                  3,    // samplestoadd
0811                                  2,    //expectedpeak
0812                                  40.,  // min energy threshold
0813                                  coef  // coefficients
0814     );
0815     result->addValues(item);
0816   }
0817   return result;
0818 }
0819 
0820 std::unique_ptr<HcalFrontEndMap> HcalHardcodeCalibrations::produceFrontEndMap(const HcalFrontEndMapRcd& rec) {
0821   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceFrontEndMap-> ...";
0822   auto const& topo = rec.get(topoTokens_[kFrontEndMap]);
0823 
0824   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0825 
0826   return dbHardcode.makeHardcodeFrontEndMap(cells);
0827 }
0828 
0829 std::unique_ptr<HcalSiPMParameters> HcalHardcodeCalibrations::produceSiPMParameters(const HcalSiPMParametersRcd& rec) {
0830   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceSiPMParameters-> ...";
0831   auto const& topo = rec.get(topoTokens_[kSiPMParameters]);
0832 
0833   auto result = std::make_unique<HcalSiPMParameters>(&topo);
0834   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0835   for (auto cell : cells) {
0836     HcalSiPMParameter item = dbHardcode.makeHardcodeSiPMParameter(cell, &topo, iLumi);
0837     result->addValues(item);
0838   }
0839   return result;
0840 }
0841 
0842 std::unique_ptr<HcalSiPMCharacteristics> HcalHardcodeCalibrations::produceSiPMCharacteristics(
0843     const HcalSiPMCharacteristicsRcd& rcd) {
0844   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceSiPMCharacteristics-> ...";
0845 
0846   return dbHardcode.makeHardcodeSiPMCharacteristics();
0847 }
0848 
0849 std::unique_ptr<HcalTPChannelParameters> HcalHardcodeCalibrations::produceTPChannelParameters(
0850     const HcalTPChannelParametersRcd& rec) {
0851   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTPChannelParameters-> ...";
0852   auto const& topo = rec.get(topoTokens_[kTPChannelParameters]);
0853 
0854   auto result = std::make_unique<HcalTPChannelParameters>(&topo);
0855   const std::vector<HcalGenericDetId>& cells = allCells(topo, dbHardcode.killHE());
0856   for (auto cell : cells) {
0857     // Thinking about Phase2 and the new FIR filter,
0858     // for now, don't put TT in TPChannelParams
0859     if (cell.subdetId() == HcalTriggerTower)
0860       continue;
0861     HcalTPChannelParameter item = dbHardcode.makeHardcodeTPChannelParameter(cell);
0862     result->addValues(item);
0863   }
0864   return result;
0865 }
0866 
0867 std::unique_ptr<HcalTPParameters> HcalHardcodeCalibrations::produceTPParameters(const HcalTPParametersRcd& rcd) {
0868   edm::LogInfo("HCAL") << "HcalHardcodeCalibrations::produceTPParameters-> ...";
0869 
0870   auto result = std::make_unique<HcalTPParameters>();
0871   dbHardcode.makeHardcodeTPParameters(*result);
0872   return result;
0873 }
0874 
0875 void HcalHardcodeCalibrations::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0876   edm::ParameterSetDescription desc;
0877   desc.add<double>("iLumi", -1.);
0878   desc.add<bool>("HBRecalibration", false);
0879   desc.add<double>("HBreCalibCutoff", 20.);
0880   desc.add<edm::FileInPath>("HBmeanenergies", edm::FileInPath("CalibCalorimetry/HcalPlugins/data/meanenergiesHB.txt"));
0881   desc.add<bool>("HERecalibration", false);
0882   desc.add<double>("HEreCalibCutoff", 20.);
0883   desc.add<edm::FileInPath>("HEmeanenergies", edm::FileInPath("CalibCalorimetry/HcalPlugins/data/meanenergiesHE.txt"));
0884   desc.add<bool>("HFRecalibration", false);
0885   desc.add<bool>("GainWidthsForTrigPrims", false);
0886   desc.add<bool>("useHBUpgrade", false);
0887   desc.add<bool>("useHEUpgrade", false);
0888   desc.add<bool>("useHFUpgrade", false);
0889   desc.add<bool>("useHOUpgrade", true);
0890   desc.add<bool>("testHFQIE10", false);
0891   desc.add<bool>("testHEPlan1", false);
0892   desc.add<bool>("killHE", false);
0893   desc.add<bool>("useLayer0Weight", false);
0894   desc.add<bool>("useIeta18depth1", true);
0895   desc.addUntracked<std::vector<std::string>>("toGet", std::vector<std::string>());
0896   desc.addUntracked<bool>("fromDDD", false);
0897 
0898   edm::ParameterSetDescription desc_hb;
0899   desc_hb.add<std::vector<double>>("gain", std::vector<double>({0.19}));
0900   desc_hb.add<std::vector<double>>("gainWidth", std::vector<double>({0.0}));
0901   desc_hb.add<double>("pedestal", 3.0);
0902   desc_hb.add<double>("pedestalWidth", 0.55);
0903   desc_hb.add<int>("zsThreshold", 8);
0904   desc_hb.add<std::vector<double>>("qieOffset", std::vector<double>({-0.49, 1.8, 7.2, 37.9}));
0905   desc_hb.add<std::vector<double>>("qieSlope", std::vector<double>({0.912, 0.917, 0.922, 0.923}));
0906   desc_hb.add<int>("qieType", 0);
0907   desc_hb.add<int>("mcShape", 125);
0908   desc_hb.add<int>("recoShape", 105);
0909   desc_hb.add<double>("photoelectronsToAnalog", 0.0);
0910   desc_hb.add<std::vector<double>>("darkCurrent", std::vector<double>({0.0}));
0911   desc_hb.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.0}));
0912   desc_hb.add<bool>("doRadiationDamage", false);
0913   desc_hb.add<double>("noiseThreshold", 0.0);
0914   desc_hb.add<double>("seedThreshold", 0.1);
0915   desc.add<edm::ParameterSetDescription>("hb", desc_hb);
0916 
0917   edm::ParameterSetDescription desc_hbRaddam;
0918   desc_hbRaddam.add<double>("temperatureBase", 20.0);
0919   desc_hbRaddam.add<double>("temperatureNew", -5.0);
0920   desc_hbRaddam.add<double>("intlumiOffset", 150);
0921   desc_hbRaddam.add<double>("depVsTemp", 0.0631);
0922   desc_hbRaddam.add<double>("intlumiToNeutrons", 3.67e8);
0923   desc_hbRaddam.add<std::vector<double>>("depVsNeutrons", {5.69e-11, 7.90e-11});
0924 
0925   edm::ParameterSetDescription desc_hbUpgrade;
0926   desc_hbUpgrade.add<std::vector<double>>("gain", std::vector<double>({0.00111111111111}));
0927   desc_hbUpgrade.add<std::vector<double>>("gainWidth", std::vector<double>({0}));
0928   desc_hbUpgrade.add<double>("pedestal", 18.0);
0929   desc_hbUpgrade.add<double>("pedestalWidth", 5.0);
0930   desc_hbUpgrade.add<int>("zsThreshold", 3);
0931   desc_hbUpgrade.add<std::vector<double>>("qieOffset", std::vector<double>({0.0, 0.0, 0.0, 0.0}));
0932   desc_hbUpgrade.add<std::vector<double>>("qieSlope", std::vector<double>({0.333, 0.333, 0.333, 0.333}));
0933   desc_hbUpgrade.add<int>("qieType", 2);
0934   desc_hbUpgrade.add<int>("mcShape", 206);
0935   desc_hbUpgrade.add<int>("recoShape", 206);
0936   desc_hbUpgrade.add<double>("photoelectronsToAnalog", 57.5);
0937   desc_hbUpgrade.add<std::vector<double>>("darkCurrent", std::vector<double>({0.055}));
0938   desc_hbUpgrade.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.26}));
0939   desc_hbUpgrade.add<bool>("doRadiationDamage", true);
0940   desc_hbUpgrade.add<edm::ParameterSetDescription>("radiationDamage", desc_hbRaddam);
0941   desc_hbUpgrade.add<double>("noiseThreshold", 0.0);
0942   desc_hbUpgrade.add<double>("seedThreshold", 0.1);
0943   desc.add<edm::ParameterSetDescription>("hbUpgrade", desc_hbUpgrade);
0944 
0945   edm::ParameterSetDescription desc_he;
0946   desc_he.add<std::vector<double>>("gain", std::vector<double>({0.23}));
0947   desc_he.add<std::vector<double>>("gainWidth", std::vector<double>({0}));
0948   desc_he.add<double>("pedestal", 3.0);
0949   desc_he.add<double>("pedestalWidth", 0.79);
0950   desc_he.add<int>("zsThreshold", 9);
0951   desc_he.add<std::vector<double>>("qieOffset", std::vector<double>({-0.38, 2.0, 7.6, 39.6}));
0952   desc_he.add<std::vector<double>>("qieSlope", std::vector<double>({0.912, 0.916, 0.92, 0.922}));
0953   desc_he.add<int>("qieType", 0);
0954   desc_he.add<int>("mcShape", 125);
0955   desc_he.add<int>("recoShape", 105);
0956   desc_he.add<double>("photoelectronsToAnalog", 0.0);
0957   desc_he.add<std::vector<double>>("darkCurrent", std::vector<double>({0.0}));
0958   desc_he.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.0}));
0959   desc_he.add<bool>("doRadiationDamage", false);
0960   desc_he.add<double>("noiseThreshold", 0.0);
0961   desc_he.add<double>("seedThreshold", 0.1);
0962   desc.add<edm::ParameterSetDescription>("he", desc_he);
0963 
0964   edm::ParameterSetDescription desc_heRaddam;
0965   desc_heRaddam.add<double>("temperatureBase", 20.0);
0966   desc_heRaddam.add<double>("temperatureNew", 5.0);
0967   desc_heRaddam.add<double>("intlumiOffset", 75);
0968   desc_heRaddam.add<double>("depVsTemp", 0.0631);
0969   desc_heRaddam.add<double>("intlumiToNeutrons", 2.92e8);
0970   desc_heRaddam.add<std::vector<double>>("depVsNeutrons", {5.69e-11, 7.90e-11});
0971 
0972   edm::ParameterSetDescription desc_heUpgrade;
0973   desc_heUpgrade.add<std::vector<double>>("gain", std::vector<double>({0.00111111111111}));
0974   desc_heUpgrade.add<std::vector<double>>("gainWidth", std::vector<double>({0}));
0975   desc_heUpgrade.add<double>("pedestal", 18.0);
0976   desc_heUpgrade.add<double>("pedestalWidth", 5.0);
0977   desc_heUpgrade.add<int>("zsThreshold", 3);
0978   desc_heUpgrade.add<std::vector<double>>("qieOffset", std::vector<double>({0.0, 0.0, 0.0, 0.0}));
0979   desc_heUpgrade.add<std::vector<double>>("qieSlope", std::vector<double>({0.333, 0.333, 0.333, 0.333}));
0980   desc_heUpgrade.add<int>("qieType", 2);
0981   desc_heUpgrade.add<int>("mcShape", 206);
0982   desc_heUpgrade.add<int>("recoShape", 206);
0983   desc_heUpgrade.add<double>("photoelectronsToAnalog", 57.5);
0984   desc_heUpgrade.add<std::vector<double>>("darkCurrent", std::vector<double>({0.055}));
0985   desc_heUpgrade.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.26}));
0986   desc_heUpgrade.add<bool>("doRadiationDamage", true);
0987   desc_heUpgrade.add<edm::ParameterSetDescription>("radiationDamage", desc_heRaddam);
0988   desc_heUpgrade.add<double>("noiseThreshold", 0.0);
0989   desc_heUpgrade.add<double>("seedThreshold", 0.1);
0990   desc.add<edm::ParameterSetDescription>("heUpgrade", desc_heUpgrade);
0991 
0992   edm::ParameterSetDescription desc_hf;
0993   desc_hf.add<std::vector<double>>("gain", std::vector<double>({0.14, 0.135}));
0994   desc_hf.add<std::vector<double>>("gainWidth", std::vector<double>({0.0, 0.0}));
0995   desc_hf.add<double>("pedestal", 3.0);
0996   desc_hf.add<double>("pedestalWidth", 0.84);
0997   desc_hf.add<int>("zsThreshold", -9999);
0998   desc_hf.add<std::vector<double>>("qieOffset", std::vector<double>({-0.87, 1.4, 7.8, -29.6}));
0999   desc_hf.add<std::vector<double>>("qieSlope", std::vector<double>({0.359, 0.358, 0.36, 0.367}));
1000   desc_hf.add<int>("qieType", 0);
1001   desc_hf.add<int>("mcShape", 301);
1002   desc_hf.add<int>("recoShape", 301);
1003   desc_hf.add<double>("photoelectronsToAnalog", 0.0);
1004   desc_hf.add<std::vector<double>>("darkCurrent", std::vector<double>({0.0}));
1005   desc_hf.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.0}));
1006   desc_hf.add<bool>("doRadiationDamage", false);
1007   desc_hf.add<double>("noiseThreshold", 0.0);
1008   desc_hf.add<double>("seedThreshold", 0.1);
1009   desc.add<edm::ParameterSetDescription>("hf", desc_hf);
1010 
1011   edm::ParameterSetDescription desc_hfUpgrade;
1012   desc_hfUpgrade.add<std::vector<double>>("gain", std::vector<double>({0.14, 0.135}));
1013   desc_hfUpgrade.add<std::vector<double>>("gainWidth", std::vector<double>({0.0, 0.0}));
1014   desc_hfUpgrade.add<double>("pedestal", 13.33);
1015   desc_hfUpgrade.add<double>("pedestalWidth", 3.33);
1016   desc_hfUpgrade.add<int>("zsThreshold", -9999);
1017   desc_hfUpgrade.add<std::vector<double>>("qieOffset", std::vector<double>({0.0697, -0.7405, 12.38, -671.9}));
1018   desc_hfUpgrade.add<std::vector<double>>("qieSlope", std::vector<double>({0.297, 0.298, 0.298, 0.313}));
1019   desc_hfUpgrade.add<int>("qieType", 1);
1020   desc_hfUpgrade.add<int>("mcShape", 301);
1021   desc_hfUpgrade.add<int>("recoShape", 301);
1022   desc_hfUpgrade.add<double>("photoelectronsToAnalog", 0.0);
1023   desc_hfUpgrade.add<std::vector<double>>("darkCurrent", std::vector<double>({0.0}));
1024   desc_hfUpgrade.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.0}));
1025   desc_hfUpgrade.add<bool>("doRadiationDamage", false);
1026   desc_hfUpgrade.add<double>("noiseThreshold", 0.0);
1027   desc_hfUpgrade.add<double>("seedThreshold", 0.1);
1028   desc.add<edm::ParameterSetDescription>("hfUpgrade", desc_hfUpgrade);
1029 
1030   edm::ParameterSetDescription desc_hfrecal;
1031   desc_hfrecal.add<std::vector<double>>("HFdepthOneParameterA", std::vector<double>());
1032   desc_hfrecal.add<std::vector<double>>("HFdepthOneParameterB", std::vector<double>());
1033   desc_hfrecal.add<std::vector<double>>("HFdepthTwoParameterA", std::vector<double>());
1034   desc_hfrecal.add<std::vector<double>>("HFdepthTwoParameterB", std::vector<double>());
1035   desc.add<edm::ParameterSetDescription>("HFRecalParameterBlock", desc_hfrecal);
1036 
1037   edm::ParameterSetDescription desc_ho;
1038   desc_ho.add<std::vector<double>>("gain", std::vector<double>({0.006, 0.0087}));
1039   desc_ho.add<std::vector<double>>("gainWidth", std::vector<double>({0.0, 0.0}));
1040   desc_ho.add<double>("pedestal", 11.0);
1041   desc_ho.add<double>("pedestalWidth", 0.57);
1042   desc_ho.add<int>("zsThreshold", 24);
1043   desc_ho.add<std::vector<double>>("qieOffset", std::vector<double>({-0.44, 1.4, 7.1, 38.5}));
1044   desc_ho.add<std::vector<double>>("qieSlope", std::vector<double>({0.907, 0.915, 0.92, 0.921}));
1045   desc_ho.add<int>("qieType", 0);
1046   desc_ho.add<int>("mcShape", 201);
1047   desc_ho.add<int>("recoShape", 201);
1048   desc_ho.add<double>("photoelectronsToAnalog", 4.0);
1049   desc_ho.add<std::vector<double>>("darkCurrent", std::vector<double>({0.0}));
1050   desc_ho.add<std::vector<double>>("noiseCorrelation", std::vector<double>({0.0}));
1051   desc_ho.add<bool>("doRadiationDamage", false);
1052   desc_ho.add<double>("noiseThreshold", 0.0);
1053   desc_ho.add<double>("seedThreshold", 0.1);
1054   desc.add<edm::ParameterSetDescription>("ho", desc_ho);
1055 
1056   edm::ParameterSetDescription validator_sipm;
1057   validator_sipm.add<int>("pixels", 1);
1058   validator_sipm.add<double>("crosstalk", 0);
1059   validator_sipm.add<double>("nonlin1", 1);
1060   validator_sipm.add<double>("nonlin2", 0);
1061   validator_sipm.add<double>("nonlin3", 0);
1062   std::vector<edm::ParameterSet> default_sipm(1);
1063   desc.addVPSet("SiPMCharacteristics", validator_sipm, default_sipm);
1064 
1065   descriptions.addDefault(desc);
1066 }