Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-16 23:24:05

0001 #include "HcalHitReconstructor.h"
0002 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0003 #include "DataFormats/Common/interface/EDCollection.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/ParameterSet/interface/EmptyGroupDescription.h"
0008 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
0009 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
0010 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
0011 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
0012 #include "CalibCalorimetry/HcalAlgos/interface/HcalDbASCIIIO.h"
0013 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0014 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0015 #include "CondFormats/DataRecord/interface/HcalFrontEndMapRcd.h"
0016 #include "CondFormats/DataRecord/interface/HcalOOTPileupCorrectionRcd.h"
0017 #include "CondFormats/DataRecord/interface/HcalOOTPileupCompatibilityRcd.h"
0018 #include "CondFormats/DataRecord/interface/HBHENegativeEFilterRcd.h"
0019 #include "CondFormats/HcalObjects/interface/HcalFrontEndMap.h"
0020 #include "CondFormats/HcalObjects/interface/OOTPileupCorrectionColl.h"
0021 #include "CondFormats/HcalObjects/interface/OOTPileupCorrData.h"
0022 #include <iostream>
0023 #include <fstream>
0024 
0025 /*  Hcal Hit reconstructor allows for CaloRecHits with status words */
0026 
0027 HcalHitReconstructor::HcalHitReconstructor(edm::ParameterSet const& conf)
0028     : reco_(conf.getParameter<bool>("correctForTimeslew"),
0029             conf.getParameter<bool>("correctForPhaseContainment"),
0030             conf.getParameter<double>("correctionPhaseNS"),
0031             consumesCollector()),
0032       det_(DetId::Hcal),
0033       inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
0034       correctTiming_(conf.getParameter<bool>("correctTiming")),
0035       setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
0036       setSaturationFlags_(conf.getParameter<bool>("setSaturationFlags")),
0037       setTimingTrustFlags_(conf.getParameter<bool>("setTimingTrustFlags")),
0038       setNegativeFlags_(conf.getParameter<bool>("setNegativeFlags")),
0039       dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
0040       firstAuxTS_(conf.getParameter<int>("firstAuxTS")),
0041       firstSample_(conf.getParameter<int>("firstSample")),
0042       samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
0043       tsFromDB_(conf.getParameter<bool>("tsFromDB")),
0044       useLeakCorrection_(conf.getParameter<bool>("useLeakCorrection")),
0045       dataOOTCorrectionName_(""),
0046       dataOOTCorrectionCategory_("Data"),
0047       mcOOTCorrectionName_(""),
0048       mcOOTCorrectionCategory_("MC"),
0049       setPileupCorrection_(nullptr) {
0050   // register for data access
0051   tok_ho_ = consumes<HODigiCollection>(inputLabel_);
0052   tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
0053   tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
0054 
0055   std::string subd = conf.getParameter<std::string>("Subdetector");
0056   //Set all FlagSetters to 0
0057   /* Important to do this!  Otherwise, if the setters are turned off,
0058      the "if (XSetter_) delete XSetter_;" commands can crash
0059   */
0060 
0061   recoParamsFromDB_ = conf.getParameter<bool>("recoParamsFromDB");
0062   //  recoParamsFromDB_ = false ; //  trun off for now.
0063 
0064   // std::cout<<"  HcalHitReconstructor   recoParamsFromDB_ "<<recoParamsFromDB_<<std::endl;
0065 
0066   hfdigibit_ = nullptr;
0067 
0068   hfS9S1_ = nullptr;
0069   hfS8S1_ = nullptr;
0070   hfPET_ = nullptr;
0071   saturationFlagSetter_ = nullptr;
0072   HFTimingTrustFlagSetter_ = nullptr;
0073   digiTimeFromDB_ = false;  // only need for HF
0074 
0075   if (setSaturationFlags_) {
0076     const edm::ParameterSet& pssat = conf.getParameter<edm::ParameterSet>("saturationParameters");
0077     saturationFlagSetter_ = new HcalADCSaturationFlag(pssat.getParameter<int>("maxADCvalue"));
0078   }
0079 
0080   if (!strcasecmp(subd.c_str(), "HO")) {
0081     subdet_ = HcalOuter;
0082     // setPileupCorrection_ = &HcalSimpleRecAlgo::setHOPileupCorrection;
0083     setPileupCorrection_ = nullptr;
0084     produces<HORecHitCollection>();
0085   } else if (!strcasecmp(subd.c_str(), "HF")) {
0086     subdet_ = HcalForward;
0087     // setPileupCorrection_ = &HcalSimpleRecAlgo::setHFPileupCorrection;
0088     setPileupCorrection_ = nullptr;
0089     digiTimeFromDB_ = conf.getParameter<bool>("digiTimeFromDB");
0090 
0091     if (setTimingTrustFlags_) {
0092       const edm::ParameterSet& pstrust = conf.getParameter<edm::ParameterSet>("hfTimingTrustParameters");
0093       HFTimingTrustFlagSetter_ = new HFTimingTrustFlag(pstrust.getParameter<int>("hfTimingTrustLevel1"),
0094                                                        pstrust.getParameter<int>("hfTimingTrustLevel2"));
0095     }
0096 
0097     if (setNoiseFlags_) {
0098       const edm::ParameterSet& psdigi = conf.getParameter<edm::ParameterSet>("digistat");
0099       const edm::ParameterSet& psTimeWin = conf.getParameter<edm::ParameterSet>("HFInWindowStat");
0100       hfdigibit_ = new HcalHFStatusBitFromDigis(psdigi, psTimeWin);
0101 
0102       const edm::ParameterSet& psS9S1 = conf.getParameter<edm::ParameterSet>("S9S1stat");
0103       hfS9S1_ = new HcalHF_S9S1algorithm(psS9S1.getParameter<std::vector<double>>("short_optimumSlope"),
0104                                          psS9S1.getParameter<std::vector<double>>("shortEnergyParams"),
0105                                          psS9S1.getParameter<std::vector<double>>("shortETParams"),
0106                                          psS9S1.getParameter<std::vector<double>>("long_optimumSlope"),
0107                                          psS9S1.getParameter<std::vector<double>>("longEnergyParams"),
0108                                          psS9S1.getParameter<std::vector<double>>("longETParams"),
0109                                          psS9S1.getParameter<int>("HcalAcceptSeverityLevel"),
0110                                          psS9S1.getParameter<bool>("isS8S1"));
0111 
0112       const edm::ParameterSet& psS8S1 = conf.getParameter<edm::ParameterSet>("S8S1stat");
0113       hfS8S1_ = new HcalHF_S9S1algorithm(psS8S1.getParameter<std::vector<double>>("short_optimumSlope"),
0114                                          psS8S1.getParameter<std::vector<double>>("shortEnergyParams"),
0115                                          psS8S1.getParameter<std::vector<double>>("shortETParams"),
0116                                          psS8S1.getParameter<std::vector<double>>("long_optimumSlope"),
0117                                          psS8S1.getParameter<std::vector<double>>("longEnergyParams"),
0118                                          psS8S1.getParameter<std::vector<double>>("longETParams"),
0119                                          psS8S1.getParameter<int>("HcalAcceptSeverityLevel"),
0120                                          psS8S1.getParameter<bool>("isS8S1"));
0121 
0122       const edm::ParameterSet& psPET = conf.getParameter<edm::ParameterSet>("PETstat");
0123       hfPET_ = new HcalHF_PETalgorithm(psPET.getParameter<std::vector<double>>("short_R"),
0124                                        psPET.getParameter<std::vector<double>>("shortEnergyParams"),
0125                                        psPET.getParameter<std::vector<double>>("shortETParams"),
0126                                        psPET.getParameter<std::vector<double>>("long_R"),
0127                                        psPET.getParameter<std::vector<double>>("longEnergyParams"),
0128                                        psPET.getParameter<std::vector<double>>("longETParams"),
0129                                        psPET.getParameter<int>("HcalAcceptSeverityLevel"),
0130                                        psPET.getParameter<std::vector<double>>("short_R_29"),
0131                                        psPET.getParameter<std::vector<double>>("long_R_29"));
0132     }
0133     produces<HFRecHitCollection>();
0134   } else if (!strcasecmp(subd.c_str(), "ZDC")) {
0135     det_ = DetId::Calo;
0136     subdet_ = HcalZDCDetId::SubdetectorId;
0137     produces<ZDCRecHitCollection>();
0138   } else if (!strcasecmp(subd.c_str(), "CALIB")) {
0139     subdet_ = HcalOther;
0140     subdetOther_ = HcalCalibration;
0141     produces<HcalCalibRecHitCollection>();
0142   } else {
0143     edm::LogWarning("Configuration") << "HcalHitReconstructor is not associated with a specific subdetector!"
0144                                      << std::endl;
0145   }
0146 
0147   // If no valid OOT pileup correction name specified,
0148   // disable the correction
0149   dataOOTCorrectionName_ = conf.getParameter<std::string>("dataOOTCorrectionName");
0150   dataOOTCorrectionCategory_ = conf.getParameter<std::string>("dataOOTCorrectionCategory");
0151   mcOOTCorrectionName_ = conf.getParameter<std::string>("mcOOTCorrectionName");
0152   mcOOTCorrectionCategory_ = conf.getParameter<std::string>("mcOOTCorrectionCategory");
0153   if (dataOOTCorrectionName_.empty() && mcOOTCorrectionName_.empty())
0154     setPileupCorrection_ = nullptr;
0155 
0156   // ES tokens
0157   htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
0158   if (tsFromDB_ || recoParamsFromDB_)
0159     paramsToken_ = esConsumes<HcalRecoParams, HcalRecoParamsRcd, edm::Transition::BeginRun>();
0160   if (digiTimeFromDB_)
0161     digiTimeToken_ = esConsumes<HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd, edm::Transition::BeginRun>();
0162   conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
0163   qualToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
0164   sevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
0165 }
0166 
0167 HcalHitReconstructor::~HcalHitReconstructor() {
0168   delete hfdigibit_;
0169 
0170   delete hfS9S1_;
0171   delete hfS8S1_;
0172   delete hfPET_;
0173   delete saturationFlagSetter_;
0174   delete HFTimingTrustFlagSetter_;
0175 }
0176 
0177 void HcalHitReconstructor::beginRun(edm::Run const& r, edm::EventSetup const& es) {
0178   const HcalTopology& htopo = es.getData(htopoToken_);
0179 
0180   if (tsFromDB_ || recoParamsFromDB_) {
0181     const HcalRecoParams& p = es.getData(paramsToken_);
0182     paramTS_ = std::make_unique<HcalRecoParams>(p);
0183     paramTS_->setTopo(&htopo);
0184 
0185     // std::cout<<" skdump in HcalHitReconstructor::beginRun   dupm RecoParams "<<std::endl;
0186     // std::ofstream skfile("skdumpRecoParamsNewFormat.txt");
0187     // HcalDbASCIIIO::dumpObject(skfile, (*paramTS_) );
0188   }
0189 
0190   if (digiTimeFromDB_) {
0191     const HcalFlagHFDigiTimeParams& p = es.getData(digiTimeToken_);
0192     hFDigiTimeParams_ = std::make_unique<HcalFlagHFDigiTimeParams>(p);
0193     hFDigiTimeParams_->setTopo(&htopo);
0194   }
0195 
0196   reco_.beginRun(es);
0197 }
0198 
0199 void HcalHitReconstructor::endRun(edm::Run const& r, edm::EventSetup const& es) { reco_.endRun(); }
0200 
0201 void HcalHitReconstructor::produce(edm::Event& e, const edm::EventSetup& eventSetup) {
0202   // get conditions
0203   const HcalDbService* conditions = &eventSetup.getData(conditionsToken_);
0204   const HcalChannelQuality* myqual = &eventSetup.getData(qualToken_);
0205   const HcalSeverityLevelComputer* mySeverity = &eventSetup.getData(sevToken_);
0206 
0207   if (useLeakCorrection_)
0208     reco_.setLeakCorrection();
0209 
0210   // GET THE BEAM CROSSING INFO HERE, WHEN WE UNDERSTAND HOW THINGS WORK.
0211   // Then, call "setBXInfo" method of the reco_ object.
0212   // Also remember to call SetBXInfo in the negative energy flag setter.
0213 
0214   if (det_ == DetId::Hcal) {
0215     //  HO ------------------------------------------------------------------
0216     if (subdet_ == HcalOuter) {
0217       edm::Handle<HODigiCollection> digi;
0218       e.getByToken(tok_ho_, digi);
0219 
0220       // create empty output
0221       auto rec = std::make_unique<HORecHitCollection>();
0222       rec->reserve(digi->size());
0223       // run the algorithm
0224       HODigiCollection::const_iterator i;
0225 
0226       // Vote on majority TS0 CapId
0227       int favorite_capid = 0;
0228       if (correctTiming_) {
0229         long capid_votes[4] = {0, 0, 0, 0};
0230         for (i = digi->begin(); i != digi->end(); i++) {
0231           capid_votes[(*i)[0].capid()]++;
0232         }
0233         for (int k = 0; k < 4; k++)
0234           if (capid_votes[k] > capid_votes[favorite_capid])
0235             favorite_capid = k;
0236       }
0237 
0238       for (i = digi->begin(); i != digi->end(); i++) {
0239         HcalDetId cell = i->id();
0240         DetId detcell = (DetId)cell;
0241         // firstSample & samplesToAdd
0242         if (tsFromDB_ || recoParamsFromDB_) {
0243           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0244           if (tsFromDB_) {
0245             firstSample_ = param_ts->firstSample();
0246             samplesToAdd_ = param_ts->samplesToAdd();
0247           }
0248           if (recoParamsFromDB_) {
0249             bool correctForTimeslew = param_ts->correctForTimeslew();
0250             bool correctForPhaseContainment = param_ts->correctForPhaseContainment();
0251             float phaseNS = param_ts->correctionPhaseNS();
0252             useLeakCorrection_ = param_ts->useLeakCorrection();
0253             correctTiming_ = param_ts->correctTiming();
0254             firstAuxTS_ = param_ts->firstAuxTS();
0255             int pileupCleaningID = param_ts->pileupCleaningID();
0256             reco_.setRecoParams(
0257                 correctForTimeslew, correctForPhaseContainment, useLeakCorrection_, pileupCleaningID, phaseNS);
0258           }
0259         }
0260 
0261         int first = firstSample_;
0262         int toadd = samplesToAdd_;
0263 
0264         if (first >= i->size() || first < 0)
0265           edm::LogWarning("Configuration")
0266               << "HcalHitReconstructor: illegal firstSample" << first << "  in subdet " << subdet_ << std::endl;
0267 
0268         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0269         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0270         if (mySeverity->dropChannel(mydigistatus->getValue()))
0271           continue;
0272         if (dropZSmarkedPassed_)
0273           if (i->zsMarkAndPass())
0274             continue;
0275 
0276         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0277         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0278         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0279         HcalCoderDb coder(*channelCoder, *shape);
0280 
0281         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0282 
0283         // Set auxiliary flag
0284         int auxflag = 0;
0285         int fTS = firstAuxTS_;
0286         if (fTS < 0)
0287           fTS = 0;  //silly protection against negative time slice values
0288         for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
0289           auxflag += (i->sample(xx).adc())
0290                      << (7 *
0291                          (xx - fTS));  // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0292         // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0293         auxflag += ((i->sample(fTS).capid()) << 28);
0294         (rec->back()).setAux(auxflag);
0295         // (rec->back()).setFlags(0);  Don't want to do this because the algorithm
0296         //                             can already set some flags
0297         // Fill Presample ADC flag
0298         if (fTS > 0)
0299           (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
0300 
0301         if (setSaturationFlags_)
0302           saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
0303         if (correctTiming_)
0304           HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
0305       }
0306       // return result
0307       e.put(std::move(rec));
0308 
0309       // HF -------------------------------------------------------------------
0310     } else if (subdet_ == HcalForward) {
0311       edm::Handle<HFDigiCollection> digi;
0312       e.getByToken(tok_hf_, digi);
0313 
0314       ///////////////////////////////////////////////////////////////// HF
0315       // create empty output
0316       auto rec = std::make_unique<HFRecHitCollection>();
0317       rec->reserve(digi->size());
0318       // run the algorithm
0319       HFDigiCollection::const_iterator i;
0320 
0321       // Vote on majority TS0 CapId
0322       int favorite_capid = 0;
0323       if (correctTiming_) {
0324         long capid_votes[4] = {0, 0, 0, 0};
0325         for (i = digi->begin(); i != digi->end(); i++) {
0326           capid_votes[(*i)[0].capid()]++;
0327         }
0328         for (int k = 0; k < 4; k++)
0329           if (capid_votes[k] > capid_votes[favorite_capid])
0330             favorite_capid = k;
0331       }
0332 
0333       for (i = digi->begin(); i != digi->end(); i++) {
0334         HcalDetId cell = i->id();
0335         DetId detcell = (DetId)cell;
0336 
0337         if (tsFromDB_ || recoParamsFromDB_) {
0338           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0339           if (tsFromDB_) {
0340             firstSample_ = param_ts->firstSample();
0341             samplesToAdd_ = param_ts->samplesToAdd();
0342           }
0343           if (recoParamsFromDB_) {
0344             bool correctForTimeslew = param_ts->correctForTimeslew();
0345             bool correctForPhaseContainment = param_ts->correctForPhaseContainment();
0346             float phaseNS = param_ts->correctionPhaseNS();
0347             useLeakCorrection_ = param_ts->useLeakCorrection();
0348             correctTiming_ = param_ts->correctTiming();
0349             firstAuxTS_ = param_ts->firstAuxTS();
0350             int pileupCleaningID = param_ts->pileupCleaningID();
0351             reco_.setRecoParams(
0352                 correctForTimeslew, correctForPhaseContainment, useLeakCorrection_, pileupCleaningID, phaseNS);
0353           }
0354         }
0355 
0356         int first = firstSample_;
0357         int toadd = samplesToAdd_;
0358 
0359         if (first >= i->size() || first < 0)
0360           edm::LogWarning("Configuration")
0361               << "HcalHitReconstructor: illegal firstSample" << first << "  in subdet " << subdet_ << std::endl;
0362 
0363         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0364         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0365         if (mySeverity->dropChannel(mydigistatus->getValue()))
0366           continue;
0367         if (dropZSmarkedPassed_)
0368           if (i->zsMarkAndPass())
0369             continue;
0370 
0371         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0372         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0373         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0374         HcalCoderDb coder(*channelCoder, *shape);
0375 
0376         // Set HFDigiTime flag values from digiTimeFromDB_
0377         if (digiTimeFromDB_ && hfdigibit_) {
0378           const HcalFlagHFDigiTimeParam* hfDTparam = hFDigiTimeParams_->getValues(detcell.rawId());
0379           hfdigibit_->resetParamsFromDB(hfDTparam->HFdigiflagFirstSample(),
0380                                         hfDTparam->HFdigiflagSamplesToAdd(),
0381                                         hfDTparam->HFdigiflagExpectedPeak(),
0382                                         hfDTparam->HFdigiflagMinEThreshold(),
0383                                         hfDTparam->HFdigiflagCoefficients());
0384         }
0385 
0386         //std::cout << "TOADDHF " << toadd << " " << first << " " << std::endl;
0387         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0388 
0389         // Set auxiliary flag
0390         int auxflag = 0;
0391         int fTS = firstAuxTS_;
0392         if (fTS < 0)
0393           fTS = 0;  // silly protection against negative time slice values
0394         for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
0395           auxflag += (i->sample(xx).adc())
0396                      << (7 *
0397                          (xx - fTS));  // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0398         // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0399         auxflag += ((i->sample(fTS).capid()) << 28);
0400         (rec->back()).setAux(auxflag);
0401 
0402         // (rec->back()).setFlags(0);  Don't want to do this because the algorithm
0403         //                             can already set some flags
0404 
0405         // Fill Presample ADC flag
0406         if (fTS > 0)
0407           (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
0408 
0409         // This calls the code for setting the HF noise bit determined from digi shape
0410         if (setNoiseFlags_)
0411           hfdigibit_->hfSetFlagFromDigi(rec->back(), *i, coder, calibrations);
0412         if (setSaturationFlags_)
0413           saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
0414         if (setTimingTrustFlags_)
0415           HFTimingTrustFlagSetter_->setHFTimingTrustFlag(rec->back(), *i);
0416         if (correctTiming_)
0417           HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
0418       }  // for (i=digi->begin(); i!=digi->end(); i++) -- loop on all HF digis
0419 
0420       // The following flags require the full set of rechits
0421       // These need to be set consecutively, so an energy check should be the first
0422       // test performed on these hits (to minimize the loop time)
0423       if (setNoiseFlags_) {
0424         // Step 1:  Set PET flag  (short fibers of |ieta|==29)
0425         // Neighbor/partner channels that are flagged by Pulse Shape algorithm (HFDigiTime)
0426         // won't be considered in these calculations
0427         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0428           int depth = i->id().depth();
0429           int ieta = i->id().ieta();
0430           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0431           if (depth == 2 || abs(ieta) == 29)
0432             hfPET_->HFSetFlagFromPET(*i, *rec, myqual, mySeverity);
0433         }
0434 
0435         // Step 2:  Set S8S1 flag (short fibers or |ieta|==29)
0436         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0437           int depth = i->id().depth();
0438           int ieta = i->id().ieta();
0439           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0440           if (depth == 2 || abs(ieta) == 29)
0441             hfS8S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
0442         }
0443 
0444         // Set 3:  Set S9S1 flag (long fibers)
0445         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0446           int depth = i->id().depth();
0447           int ieta = i->id().ieta();
0448           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0449           if (depth == 1 && abs(ieta) != 29)
0450             hfS9S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
0451         }
0452       }
0453 
0454       // return result
0455       e.put(std::move(rec));
0456     } else if (subdet_ == HcalOther && subdetOther_ == HcalCalibration) {
0457       edm::Handle<HcalCalibDigiCollection> digi;
0458       e.getByToken(tok_calib_, digi);
0459 
0460       // create empty output
0461       auto rec = std::make_unique<HcalCalibRecHitCollection>();
0462       rec->reserve(digi->size());
0463       // run the algorithm
0464       int first = firstSample_;
0465       int toadd = samplesToAdd_;
0466 
0467       HcalCalibDigiCollection::const_iterator i;
0468       for (i = digi->begin(); i != digi->end(); i++) {
0469         HcalCalibDetId cell = i->id();
0470         //  HcalDetId cellh = i->id();
0471         DetId detcell = (DetId)cell;
0472         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0473         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0474         if (mySeverity->dropChannel(mydigistatus->getValue()))
0475           continue;
0476         if (dropZSmarkedPassed_)
0477           if (i->zsMarkAndPass())
0478             continue;
0479 
0480         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0481         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0482         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0483         HcalCoderDb coder(*channelCoder, *shape);
0484 
0485         // firstSample & samplesToAdd
0486         if (tsFromDB_) {
0487           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0488           first = param_ts->firstSample();
0489           toadd = param_ts->samplesToAdd();
0490         }
0491         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0492 
0493         /*
0494       // Flag setting not available for calibration rechits
0495     // Set auxiliary flag
0496     int auxflag=0;
0497         int fTS = firstAuxTS_;
0498     for (int xx=fTS; xx<fTS+4 && xx<i->size();++xx)
0499       auxflag+=(i->sample(xx).adc())<<(7*(xx-fTS)); // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0500     // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0501     auxflag+=((i->sample(fTS).capid())<<28);
0502     (rec->back()).setAux(auxflag);
0503 
0504     (rec->back()).setFlags(0); // Not yet implemented for HcalCalibRecHit
0505     */
0506       }
0507       // return result
0508       e.put(std::move(rec));
0509     }
0510   }
0511   //DL  delete myqual;
0512 }  // void HcalHitReconstructor::produce(...)
0513 
0514 void HcalHitReconstructor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0515   edm::ParameterSetDescription desc;
0516   desc.add<bool>("correctForTimeslew", false);
0517   desc.add<bool>("correctForPhaseContainment", false);
0518   desc.add<double>("correctionPhaseNS", 13.0);
0519   desc.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
0520   desc.add<bool>("correctTiming", true);
0521   desc.add<bool>("dropZSmarkedPassed", true);
0522   desc.add<int>("firstAuxTS", 1);
0523   desc.add<int>("firstSample", 2);
0524   desc.add<int>("samplesToAdd", 1);
0525   desc.add<bool>("tsFromDB", true);
0526   desc.add<bool>("useLeakCorrection", false);
0527   desc.add<bool>("recoParamsFromDB", true);
0528   desc.add<bool>("setNegativeFlags", false);
0529 
0530   edm::ParameterSetDescription saturationParametersDesc;
0531   saturationParametersDesc.add<int>("maxADCvalue", 127);
0532   desc.add<edm::ParameterSetDescription>("saturationParameters", saturationParametersDesc);
0533 
0534   desc.add<bool>("setSaturationFlags", true);
0535   desc.add<std::string>("Subdetector", "HF");
0536   desc.add<bool>("digiTimeFromDB", false);
0537 
0538   edm::ParameterSetDescription hfTimingTrustParametersDesc;
0539   hfTimingTrustParametersDesc.add<int>("hfTimingTrustLevel1", 1);
0540   hfTimingTrustParametersDesc.add<int>("hfTimingTrustLevel2", 4);
0541   desc.add<edm::ParameterSetDescription>("hfTimingTrustParameters", hfTimingTrustParametersDesc);
0542 
0543   desc.add<bool>("setTimingTrustFlags", true);
0544   desc.add<bool>("setNoiseFlags", true);
0545 
0546   edm::ParameterSetDescription digiStatDesc;
0547   HcalHFStatusBitFromDigis::fillHFDigiTimeParamsDesc(digiStatDesc);
0548   desc.add<edm::ParameterSetDescription>("digistat", digiStatDesc);
0549 
0550   edm::ParameterSetDescription hfInWindowStatDesc;
0551   HcalHFStatusBitFromDigis::fillHFTimeInWindowParamsDesc(hfInWindowStatDesc);
0552   desc.add<edm::ParameterSetDescription>("HFInWindowStat", hfInWindowStatDesc);
0553 
0554   {
0555     edm::ParameterSetDescription s9s1StatDesc;
0556     s9s1StatDesc.add<std::vector<double>>("short_optimumSlope",
0557                                           {-99999,
0558                                            0.0164905,
0559                                            0.0238698,
0560                                            0.0321383,
0561                                            0.041296,
0562                                            0.0513428,
0563                                            0.0622789,
0564                                            0.0741041,
0565                                            0.0868186,
0566                                            0.100422,
0567                                            0.135313,
0568                                            0.136289,
0569                                            0.0589927});
0570     s9s1StatDesc.add<std::vector<double>>(
0571         "shortEnergyParams",
0572         {35.1773, 35.37, 35.7933, 36.4472, 37.3317, 38.4468, 39.7925, 41.3688, 43.1757, 45.2132, 47.4813, 49.98, 52.7093});
0573     s9s1StatDesc.add<std::vector<double>>("shortETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0574     s9s1StatDesc.add<std::vector<double>>("long_optimumSlope",
0575                                           {-99999,
0576                                            0.0164905,
0577                                            0.0238698,
0578                                            0.0321383,
0579                                            0.041296,
0580                                            0.0513428,
0581                                            0.0622789,
0582                                            0.0741041,
0583                                            0.0868186,
0584                                            0.100422,
0585                                            0.135313,
0586                                            0.136289,
0587                                            0.0589927});
0588     s9s1StatDesc.add<std::vector<double>>(
0589         "longEnergyParams", {43.5, 45.7, 48.32, 51.36, 54.82, 58.7, 63.0, 67.72, 72.86, 78.42, 84.4, 90.8, 97.62});
0590     s9s1StatDesc.add<std::vector<double>>("longETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0591     s9s1StatDesc.add<int>("HcalAcceptSeverityLevel", 9);
0592     s9s1StatDesc.add<bool>("isS8S1", false);
0593     desc.add<edm::ParameterSetDescription>("S9S1stat", s9s1StatDesc);
0594   }
0595 
0596   {
0597     edm::ParameterSetDescription s8s1StatDesc;
0598     s8s1StatDesc.add<std::vector<double>>(
0599         "short_optimumSlope", {0.30, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10});
0600     s8s1StatDesc.add<std::vector<double>>("shortEnergyParams",
0601                                           {40, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100});
0602     s8s1StatDesc.add<std::vector<double>>("shortETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0603     s8s1StatDesc.add<std::vector<double>>(
0604         "long_optimumSlope", {0.30, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10, 0.10});
0605     s8s1StatDesc.add<std::vector<double>>("longEnergyParams",
0606                                           {40, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100, 100});
0607     s8s1StatDesc.add<std::vector<double>>("longETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0608     s8s1StatDesc.add<int>("HcalAcceptSeverityLevel", 9);
0609     s8s1StatDesc.add<bool>("isS8S1", true);
0610     desc.add<edm::ParameterSetDescription>("S8S1stat", s8s1StatDesc);
0611   }
0612 
0613   {
0614     edm::ParameterSetDescription petStatDesc;
0615     petStatDesc.add<std::vector<double>>("short_R", {0.8});
0616     petStatDesc.add<std::vector<double>>(
0617         "shortEnergyParams",
0618         {35.1773, 35.37, 35.7933, 36.4472, 37.3317, 38.4468, 39.7925, 41.3688, 43.1757, 45.2132, 47.4813, 49.98, 52.7093});
0619     petStatDesc.add<std::vector<double>>("shortETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0620     petStatDesc.add<std::vector<double>>("long_R", {0.98});
0621     petStatDesc.add<std::vector<double>>(
0622         "longEnergyParams", {43.5, 45.7, 48.32, 51.36, 54.82, 58.7, 63.0, 67.72, 72.86, 78.42, 84.4, 90.8, 97.62});
0623     petStatDesc.add<std::vector<double>>("longETParams", {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0});
0624     petStatDesc.add<std::vector<double>>("short_R_29", {0.8});
0625     petStatDesc.add<std::vector<double>>("long_R_29", {0.8});
0626     petStatDesc.add<int>("HcalAcceptSeverityLevel", 9);
0627     desc.add<edm::ParameterSetDescription>("PETstat", petStatDesc);
0628   }
0629 
0630   desc.add<std::string>("dataOOTCorrectionName", std::string(""));
0631   desc.add<std::string>("dataOOTCorrectionCategory", "Data");
0632   desc.add<std::string>("mcOOTCorrectionName", "");
0633   desc.add<std::string>("mcOOTCorrectionCategory", "MC");
0634   descriptions.addWithDefaultLabel(desc);
0635 }